Skip to content

Commit 95cfdd2

Browse files
authored
Merge pull request #820 from jhdark/cyl_and_sph_avg_quantities
Average surface and volume derived quantities in Cylindrical and Spherical
2 parents 8d8f7d9 + 10fc3b0 commit 95cfdd2

File tree

6 files changed

+439
-6
lines changed

6 files changed

+439
-6
lines changed

festim/__init__.py

Lines changed: 10 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -74,14 +74,22 @@
7474
)
7575
from .exports.derived_quantities.hydrogen_flux import HydrogenFlux
7676
from .exports.derived_quantities.thermal_flux import ThermalFlux
77-
from .exports.derived_quantities.average_volume import AverageVolume
77+
from .exports.derived_quantities.average_volume import (
78+
AverageVolume,
79+
AverageVolumeCylindrical,
80+
AverageVolumeSpherical,
81+
)
7882
from .exports.derived_quantities.maximum_volume import MaximumVolume
7983
from .exports.derived_quantities.minimum_volume import MinimumVolume
8084
from .exports.derived_quantities.minimum_surface import MinimumSurface
8185
from .exports.derived_quantities.maximum_surface import MaximumSurface
8286
from .exports.derived_quantities.total_surface import TotalSurface
8387
from .exports.derived_quantities.total_volume import TotalVolume
84-
from .exports.derived_quantities.average_surface import AverageSurface
88+
from .exports.derived_quantities.average_surface import (
89+
AverageSurface,
90+
AverageSurfaceCylindrical,
91+
AverageSurfaceSpherical,
92+
)
8593
from .exports.derived_quantities.point_value import PointValue
8694
from .exports.derived_quantities.adsorbed_hydrogen import AdsorbedHydrogen
8795

festim/exports/derived_quantities/average_surface.py

Lines changed: 57 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -30,7 +30,7 @@ def __init__(self, field, surface) -> None:
3030

3131
@property
3232
def allowed_meshes(self):
33-
return ["cartesian"]
33+
return ["cartesian", "spherical"]
3434

3535
@property
3636
def export_unit(self):
@@ -51,3 +51,59 @@ def compute(self):
5151
return f.assemble(self.function * self.ds(self.surface)) / f.assemble(
5252
1 * self.ds(self.surface)
5353
)
54+
55+
56+
class AverageSurfaceCylindrical(AverageSurface):
57+
"""
58+
Computes the average value of a field on a given surface
59+
int(f ds) / int (1 * ds)
60+
ds is the surface measure in cylindrical coordinates.
61+
ds = r dr dtheta
62+
63+
Args:
64+
field (str, int): the field ("solute", 0, 1, "T", "retention")
65+
surface (int): the surface id
66+
67+
Attributes:
68+
field (str, int): the field ("solute", 0, 1, "T", "retention")
69+
surface (int): the surface id
70+
title (str): the title of the derived quantity
71+
show_units (bool): show the units in the title in the derived quantities
72+
file
73+
function (dolfin.function.function.Function): the solution function of
74+
the field
75+
r (ufl.indexed.Indexed): the radius of the cylinder
76+
77+
Notes:
78+
Units are in H/m3 for hydrogen concentration and K for temperature
79+
"""
80+
81+
def __init__(self, field, surface) -> None:
82+
super().__init__(field=field, surface=surface)
83+
self.r = None
84+
85+
@property
86+
def allowed_meshes(self):
87+
return ["cylindrical"]
88+
89+
def compute(self):
90+
91+
if self.r is None:
92+
mesh = (
93+
self.function.function_space().mesh()
94+
) # get the mesh from the function
95+
rthetaz = f.SpatialCoordinate(mesh) # get the coordinates from the mesh
96+
self.r = rthetaz[0] # only care about r here
97+
98+
# dS_z = r dr dtheta , assuming axisymmetry dS_z = theta r dr
99+
# dS_r = r dz dtheta , assuming axisymmetry dS_r = theta r dz
100+
# in both cases the expression with self.dx is the same
101+
102+
avg_surf = f.assemble(
103+
self.function * self.r * self.ds(self.surface)
104+
) / f.assemble(1 * self.r * self.ds(self.surface))
105+
106+
return avg_surf
107+
108+
109+
AverageSurfaceSpherical = AverageSurface

festim/exports/derived_quantities/average_volume.py

Lines changed: 92 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -19,6 +19,7 @@ class AverageVolume(VolumeQuantity):
1919
file
2020
function (dolfin.function.function.Function): the solution function of
2121
the field
22+
r (ufl.indexed.Indexed): the radius of the cylinder
2223
2324
.. note::
2425
Units are in H/m3 for hydrogen concentration and K for temperature
@@ -50,3 +51,94 @@ def compute(self):
5051
return f.assemble(self.function * self.dx(self.volume)) / f.assemble(
5152
1 * self.dx(self.volume)
5253
)
54+
55+
56+
class AverageVolumeCylindrical(AverageVolume):
57+
"""
58+
Computes the average value of a field in a given volume
59+
int(f dx) / int (1 * dx)
60+
dx is the volume measure in cylindrical coordinates.
61+
dx = r dr dz dtheta
62+
63+
Note: for particle fluxes J is given in H/s, for heat fluxes J is given in W
64+
65+
Args:
66+
field (str, int): the field ("solute", 0, 1, "T", "retention")
67+
volume (int): the volume id
68+
69+
Attributes:
70+
field (str, int): the field ("solute", 0, 1, "T", "retention")
71+
volume (int): the volume id
72+
title (str): the title of the derived quantity
73+
show_units (bool): show the units in the title in the derived quantities
74+
file
75+
function (dolfin.function.function.Function): the solution function of
76+
the field
77+
r (ufl.indexed.Indexed): the radius of the sphere
78+
"""
79+
80+
def __init__(self, field, volume) -> None:
81+
super().__init__(field=field, volume=volume)
82+
self.r = None
83+
84+
@property
85+
def allowed_meshes(self):
86+
return ["cylindrical"]
87+
88+
def compute(self):
89+
90+
if self.r is None:
91+
mesh = (
92+
self.function.function_space().mesh()
93+
) # get the mesh from the function
94+
rthetaz = f.SpatialCoordinate(mesh) # get the coordinates from the mesh
95+
self.r = rthetaz[0] # only care about r here
96+
97+
avg_vol = f.assemble(
98+
self.function * self.r * self.dx(self.volume)
99+
) / f.assemble(1 * self.r * self.dx(self.volume))
100+
101+
return avg_vol
102+
103+
104+
class AverageVolumeSpherical(AverageVolume):
105+
"""
106+
Computes the average value of a field in a given volume
107+
int(f dx) / int (1 * dx)
108+
dx is the volume measure in cylindrical coordinates.
109+
dx = rho dtheta dphi
110+
111+
Note: for particle fluxes J is given in H/s, for heat fluxes J is given in W
112+
113+
Args:
114+
field (str, int): the field ("solute", 0, 1, "T", "retention")
115+
volume (int): the volume id
116+
title (str): the title of the derived quantity
117+
show_units (bool): show the units in the title in the derived quantities
118+
file
119+
function (dolfin.function.function.Function): the solution function of
120+
the field
121+
"""
122+
123+
def __init__(self, field, volume) -> None:
124+
super().__init__(field=field, volume=volume)
125+
self.r = None
126+
127+
@property
128+
def allowed_meshes(self):
129+
return ["spherical"]
130+
131+
def compute(self):
132+
133+
if self.r is None:
134+
mesh = (
135+
self.function.function_space().mesh()
136+
) # get the mesh from the function
137+
rthetaphi = f.SpatialCoordinate(mesh) # get the coordinates from the mesh
138+
self.r = rthetaphi[0] # only care about r here
139+
140+
avg_vol = f.assemble(
141+
self.function * self.r**2 * self.dx(self.volume)
142+
) / f.assemble(1 * self.r**2 * self.dx(self.volume))
143+
144+
return avg_vol

test/simulation/test_initialise.py

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -117,7 +117,6 @@ def test_TXTExport_times_added_to_milestones(tmpdir):
117117
F.SurfaceFlux(field="solute", surface=1),
118118
F.TotalVolume(field="solute", volume=1),
119119
F.TotalSurface(field="solute", surface=1),
120-
F.AverageSurface(field="solute", surface=1),
121120
F.AverageVolume(field="solute", volume=1),
122121
F.HydrogenFlux(surface=1),
123122
F.ThermalFlux(surface=1),

test/unit/test_exports/test_derived_quantities/test_average_surface.py

Lines changed: 153 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,14 @@
1-
from festim import AverageSurface
1+
from festim import (
2+
AverageSurface,
3+
AverageSurfaceCylindrical,
4+
AverageSurfaceSpherical,
5+
x,
6+
y,
7+
)
28
import fenics as f
39
import pytest
10+
import numpy as np
11+
from sympy.printing import ccode
412

513

614
@pytest.mark.parametrize("field, surface", [("solute", 1), ("T", 2)])
@@ -46,3 +54,147 @@ def test_h_average(self):
4654
)
4755
computed = self.my_average.compute()
4856
assert computed == expected
57+
58+
59+
@pytest.mark.parametrize("radius", [1, 4])
60+
@pytest.mark.parametrize("r0", [3, 5])
61+
@pytest.mark.parametrize("height", [2, 7])
62+
def test_compute_cylindrical(r0, radius, height):
63+
"""
64+
Test that AverageSurfaceCylindrical computes the value correctly on a hollow cylinder
65+
66+
Args:
67+
r0 (float): internal radius
68+
radius (float): cylinder radius
69+
height (float): cylinder height
70+
"""
71+
# creating a mesh with FEniCS
72+
r1 = r0 + radius
73+
z0, z1 = 0, height
74+
75+
mesh_fenics = f.RectangleMesh(f.Point(r0, z0), f.Point(r1, z1), 10, 10)
76+
77+
top_surface = f.CompiledSubDomain(
78+
f"on_boundary && near(x[1], {z1}, tol)", tol=1e-14
79+
)
80+
surface_markers = f.MeshFunction(
81+
"size_t", mesh_fenics, mesh_fenics.topology().dim() - 1
82+
)
83+
surface_markers.set_all(0)
84+
ds = f.Measure("ds", domain=mesh_fenics, subdomain_data=surface_markers)
85+
# Surface ids
86+
top_id = 2
87+
top_surface.mark(surface_markers, top_id)
88+
89+
my_export = AverageSurfaceCylindrical("solute", top_id)
90+
V = f.FunctionSpace(mesh_fenics, "P", 1)
91+
c_fun = lambda r: 2 * r
92+
93+
expr = f.Expression(
94+
ccode(c_fun(x)),
95+
degree=1,
96+
)
97+
my_export.function = f.interpolate(expr, V)
98+
my_export.ds = ds
99+
100+
expected_value = 4 / 3 * (r1**3 - r0**3) / (r1**2 - r0**2)
101+
102+
computed_value = my_export.compute()
103+
104+
assert np.isclose(computed_value, expected_value)
105+
106+
107+
@pytest.mark.parametrize("radius", [2, 4])
108+
@pytest.mark.parametrize("r0", [3, 5])
109+
def test_compute_spherical(r0, radius):
110+
"""
111+
Test that AverageSurfaceSpherical computes the average value correctly
112+
on a hollow sphere
113+
114+
Args:
115+
r0 (float): internal radius
116+
radius (float): sphere radius
117+
"""
118+
# creating a mesh with FEniCS
119+
r1 = r0 + radius
120+
mesh_fenics = f.IntervalMesh(10, r0, r1)
121+
122+
# marking physical groups (volumes and surfaces)
123+
right_surface = f.CompiledSubDomain(
124+
f"on_boundary && near(x[0], {r1}, tol)", tol=1e-14
125+
)
126+
surface_markers = f.MeshFunction(
127+
"size_t", mesh_fenics, mesh_fenics.topology().dim() - 1
128+
)
129+
surface_markers.set_all(0)
130+
# Surface ids
131+
right_id = 2
132+
right_surface.mark(surface_markers, right_id)
133+
ds = f.Measure("ds", domain=mesh_fenics, subdomain_data=surface_markers)
134+
135+
my_export = AverageSurfaceSpherical("solute", right_id)
136+
V = f.FunctionSpace(mesh_fenics, "P", 1)
137+
c_fun = lambda r: 4 * r
138+
expr = f.Expression(
139+
ccode(c_fun(x)),
140+
degree=1,
141+
)
142+
my_export.function = f.interpolate(expr, V)
143+
144+
my_export.ds = ds
145+
146+
expected_value = 4 * r1
147+
148+
computed_value = my_export.compute()
149+
150+
assert np.isclose(computed_value, expected_value)
151+
152+
153+
def test_average_surface_cylindrical_title_no_units_solute():
154+
"""A simple test to check that the title is set correctly in
155+
festim.AverageSurfaceCylindrical with a solute field without units"""
156+
157+
my_export = AverageSurfaceCylindrical("solute", 4)
158+
assert my_export.title == "Average solute surface 4 (H m-3)"
159+
160+
161+
def test_average_surface_cylindrical_title_no_units_temperature():
162+
"""A simple test to check that the title is set correctly in
163+
festim.AverageSurfaceCylindrical with a T field without units"""
164+
165+
my_export = AverageSurfaceCylindrical("T", 5)
166+
assert my_export.title == "Average T surface 5 (K)"
167+
168+
169+
def test_average_surface_spherical_title_no_units_solute():
170+
"""A simple test to check that the title is set correctly in
171+
festim.AverageSurfaceSpherical with a solute field without units"""
172+
173+
my_export = AverageSurfaceSpherical("solute", 6)
174+
assert my_export.title == "Average solute surface 6 (H m-3)"
175+
176+
177+
def test_average_surface_spherical_title_no_units_temperature():
178+
"""A simple test to check that the title is set correctly in
179+
festim.AverageSurfaceSpherical with a T field without units"""
180+
181+
my_export = AverageSurfaceSpherical("T", 9)
182+
assert my_export.title == "Average T surface 9 (K)"
183+
184+
185+
def test_avg_surf_cylindrical_allow_meshes():
186+
"""A simple test to check cylindrical meshes are the only
187+
meshes allowed when using AverageSurfaceCylindrical"""
188+
189+
my_export = AverageSurfaceCylindrical("solute", 2)
190+
191+
assert my_export.allowed_meshes == ["cylindrical"]
192+
193+
194+
def test_avg_surf_spherical_allow_meshes():
195+
"""A simple test to check spherical meshes are one of the
196+
meshes allowed when using AverageSurfaceSpherical"""
197+
198+
my_export = AverageSurfaceSpherical("T", 6)
199+
200+
assert "spherical" in my_export.allowed_meshes

0 commit comments

Comments
 (0)