Skip to content

Commit

Permalink
Merge pull request #689 from RemDelaporteMathurin/cylindrical_flux
Browse files Browse the repository at this point in the history
Cylindrical and spherical versions of SurfaceFlux
  • Loading branch information
RemDelaporteMathurin authored Mar 25, 2024
2 parents dde42a7 + 7d841f0 commit 583572d
Show file tree
Hide file tree
Showing 6 changed files with 363 additions and 21 deletions.
8 changes: 8 additions & 0 deletions docs/source/api/exports.rst
Original file line number Diff line number Diff line change
Expand Up @@ -47,6 +47,14 @@ Exports
:members:
:show-inheritance:

.. autoclass:: SurfaceFluxCylindrical
:members:
:show-inheritance:

.. autoclass:: SurfaceFluxSpherical
:members:
:show-inheritance:

.. autoclass:: HydrogenFlux
:members:
:show-inheritance:
Expand Down
6 changes: 5 additions & 1 deletion festim/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -67,7 +67,11 @@
VolumeQuantity,
SurfaceQuantity,
)
from .exports.derived_quantities.surface_flux import SurfaceFlux
from .exports.derived_quantities.surface_flux import (
SurfaceFlux,
SurfaceFluxCylindrical,
SurfaceFluxSpherical,
)
from .exports.derived_quantities.hydrogen_flux import HydrogenFlux
from .exports.derived_quantities.thermal_flux import ThermalFlux
from .exports.derived_quantities.average_volume import AverageVolume
Expand Down
146 changes: 144 additions & 2 deletions festim/exports/derived_quantities/surface_flux.py
Original file line number Diff line number Diff line change
@@ -1,5 +1,6 @@
from festim import SurfaceQuantity, k_B
import fenics as f
import numpy as np


class SurfaceFlux(SurfaceQuantity):
Expand All @@ -22,14 +23,17 @@ def __init__(self, field, surface) -> None:
super().__init__(field=field, surface=surface)
self.title = "Flux surface {}: {}".format(self.surface, self.field)

def compute(self, soret=False):
@property
def prop(self):
field_to_prop = {
"0": self.D,
"solute": self.D,
0: self.D,
"T": self.thermal_cond,
}
self.prop = field_to_prop[self.field]
return field_to_prop[self.field]

def compute(self, soret=False):
flux = f.assemble(
self.prop * f.dot(f.grad(self.function), self.n) * self.ds(self.surface)
)
Expand All @@ -43,3 +47,141 @@ def compute(self, soret=False):
* self.ds(self.surface)
)
return flux


class SurfaceFluxCylindrical(SurfaceFlux):
"""
Object to compute the flux J of a field u through a surface
J = integral(-prop * grad(u) . n ds)
where prop is the property of the field (D, thermal conductivity, etc)
u is the field
n is the normal vector of the surface
ds is the surface measure in cylindrical coordinates.
ds = r dr dtheta or ds = r dz dtheta
Note: for particle fluxes J is given in H/s, for heat fluxes J is given in W
Args:
field (str, int): the field ("solute", 0, 1, "T", "retention")
surface (int): the surface id
azimuth_range (tuple, optional): Range of the azimuthal angle
(theta) needs to be between 0 and 2 pi. Defaults to (0, 2 * np.pi).
"""

def __init__(self, field, surface, azimuth_range=(0, 2 * np.pi)) -> None:
super().__init__(field, surface)
self.r = None
self.azimuth_range = azimuth_range

@property
def azimuth_range(self):
return self._azimuth_range

@azimuth_range.setter
def azimuth_range(self, value):
if value[0] < 0 or value[1] > 2 * np.pi:
raise ValueError("Azimuthal range must be between 0 and pi")
self._azimuth_range = value

def compute(self, soret=False):
if soret:
raise NotImplementedError(
"Soret effect not implemented for cylindrical coordinates"
)

if self.r is None:
mesh = (
self.function.function_space().mesh()
) # get the mesh from the function
rthetaz = f.SpatialCoordinate(mesh) # get the coordinates from the mesh
self.r = rthetaz[0] # only care about r here

# dS_z = r dr dtheta , assuming axisymmetry dS_z = theta r dr
# dS_r = r dz dtheta , assuming axisymmetry dS_r = theta r dz
# in both cases the expression with self.ds is the same
# we assume full cylinder theta = 2 pi
flux = f.assemble(
self.prop
* self.r
* f.dot(f.grad(self.function), self.n)
* self.ds(self.surface)
)
flux *= self.azimuth_range[1] - self.azimuth_range[0]
return flux


class SurfaceFluxSpherical(SurfaceFlux):
"""
Object to compute the flux J of a field u through a surface
J = integral(-prop * grad(u) . n ds)
where prop is the property of the field (D, thermal conductivity, etc)
u is the field
n is the normal vector of the surface
ds is the surface measure in spherical coordinates.
ds = r^2 sin(theta) dtheta dphi
Note: for particle fluxes J is given in H/s, for heat fluxes J is given in W
Args:
field (str, int): the field ("solute", 0, 1, "T", "retention")
surface (int): the surface id
azimuth_range (tuple, optional): Range of the azimuthal angle
(phi) needs to be between 0 and pi. Defaults to (0, np.pi).
polar_range (tuple, optional): Range of the polar angle
(theta) needs to be between - pi and pi. Defaults to (-np.pi, np.pi).
"""

def __init__(
self, field, surface, azimuth_range=(0, np.pi), polar_range=(-np.pi, np.pi)
) -> None:
super().__init__(field, surface)
self.r = None
self.polar_range = polar_range
self.azimuth_range = azimuth_range

@property
def polar_range(self):
return self._polar_range

@polar_range.setter
def polar_range(self, value):
if value[0] < -np.pi or value[1] > np.pi:
raise ValueError("Polar range must be between - pi and pi")
self._polar_range = value

@property
def azimuth_range(self):
return self._azimuth_range

@azimuth_range.setter
def azimuth_range(self, value):
if value[0] < 0 or value[1] > np.pi:
raise ValueError("Azimuthal range must be between 0 and pi")
self._azimuth_range = value

def compute(self, soret=False):
if soret:
raise NotImplementedError(
"Soret effect not implemented for spherical coordinates"
)

if self.r is None:
mesh = (
self.function.function_space().mesh()
) # get the mesh from the function
rthetaphi = f.SpatialCoordinate(mesh) # get the coordinates from the mesh
self.r = rthetaphi[0] # only care about r here

# dS_r = r^2 sin(theta) dtheta dphi
# integral(f dS_r) = integral(f r^2 sin(theta) dtheta dphi)
# = (phi2 - phi1) * (-cos(theta2) + cos(theta1)) * f r^2
flux = f.assemble(
self.prop
* self.r**2
* f.dot(f.grad(self.function), self.n)
* self.ds(self.surface)
)
flux *= (self.polar_range[1] - self.polar_range[0]) * (
-np.cos(self.azimuth_range[1]) + np.cos(self.azimuth_range[0])
)
return flux
39 changes: 23 additions & 16 deletions festim/generic_simulation.py
Original file line number Diff line number Diff line change
Expand Up @@ -297,25 +297,32 @@ def initialise(self):

self.h_transport_problem.initialise(self.mesh, self.materials, self.dt)

# raise warning if SurfaceFlux is used with non-cartesian meshes
# raise warning if the derived quantities don't match the type of mesh
# eg. SurfaceFlux is used with cylindrical mesh
all_types_quantities = [
festim.MaximumSurface,
festim.MinimumSurface,
festim.MaximumVolume,
festim.MinimumVolume,
festim.PointValue,
] # these quantities can be used with any mesh
allowed_quantities = {
"cartesian": [
festim.SurfaceFlux,
festim.AverageSurface,
festim.AverageVolume,
]
+ all_types_quantities,
"cylindrical": [festim.SurfaceFluxCylindrical] + all_types_quantities,
"spherical": [festim.SurfaceFluxSpherical] + all_types_quantities,
}

for export in self.exports:
if isinstance(export, festim.DerivedQuantities):
allowed_quantities = (
festim.MaximumSurface,
festim.MinimumSurface,
festim.MaximumVolume,
festim.MinimumVolume,
festim.PointValue,
)
if any(
[
not isinstance(q, allowed_quantities)
for q in export.derived_quantities
]
):
if self.mesh.type != "cartesian":
for q in export.derived_quantities:
if not isinstance(q, tuple(allowed_quantities[self.mesh.type])):
warnings.warn(
"Some derived quantities may not work as intended for non-cartesian meshes"
f"{type(q)} may not work as intended for {self.mesh.type} meshes"
)
self.exports.initialise_derived_quantities(
self.mesh.dx, self.mesh.ds, self.materials
Expand Down
2 changes: 1 addition & 1 deletion test/simulation/test_initialise.py
Original file line number Diff line number Diff line change
Expand Up @@ -137,5 +137,5 @@ def test_cartesian_and_surface_flux_warning(quantity, sys):
my_model.exports = [derived_quantities]

# test
with pytest.warns(UserWarning, match="Some derived quantities .* non-cartesian"):
with pytest.warns(UserWarning, match=f"may not work as intended for {sys} meshes"):
my_model.initialise()
Loading

0 comments on commit 583572d

Please sign in to comment.