Skip to content

Commit

Permalink
Added 'local_rhs' option to 'local_project' and 'calculate_and_print_…
Browse files Browse the repository at this point in the history
…flow_properties' functions as recommended by Jørgen
  • Loading branch information
johannesring committed Sep 7, 2023
1 parent a28b0aa commit 6cfcc39
Showing 1 changed file with 11 additions and 4 deletions.
15 changes: 11 additions & 4 deletions src/fsipy/simulations/simulation_common.py
Original file line number Diff line number Diff line change
Expand Up @@ -211,13 +211,15 @@ def peval(f: Function, x: Union[float, np.ndarray]) -> np.ndarray:
return yglob


def local_project(f: Function, V: FunctionSpace) -> Function:
def local_project(f: Function, V: FunctionSpace, local_rhs: bool = False) -> Function:
"""
Project a given function 'f' onto a finite element function space 'V' in a reusable way.
Args:
f (Function): The function to be projected onto 'V'.
V (FunctionSpace): The finite element function space.
local_rhs (bool, optional): If True, solve using a local right-hand side assembly.
If False (default), solve using a global right-hand side assembly.
Returns:
Function: The projected solution in 'V'.
Expand All @@ -229,13 +231,16 @@ def local_project(f: Function, V: FunctionSpace) -> Function:
solver = LocalSolver(a_proj, b_proj)
solver.factorize()
projected_u = Function(V)
solver.solve_local_rhs(projected_u)
if local_rhs:
solver.solve_local_rhs(projected_u)
else:
solver.solve_global_rhs(projected_u)

return projected_u


def calculate_and_print_flow_properties(dt: float, mesh: Mesh, v: Function, inlet_area: float, mu_f: float,
n: Expression, dsi: Measure) -> None:
n: Expression, dsi: Measure, local_rhs: bool = False) -> None:
"""
Calculate and print flow properties.
Expand All @@ -247,13 +252,15 @@ def calculate_and_print_flow_properties(dt: float, mesh: Mesh, v: Function, inle
mu_f (float): Fluid dynamic viscosity.
n (dolfin.Expression): FacetNormal expression.
dsi (dolfin.Measure): Measure for inlet boundary.
local_rhs (bool, optional): If True, solve using a local right-hand side assembly.
If False (default), solve using a global right-hand side assembly.
Returns:
None
"""
# Calculate the DG vector of velocity magnitudes
DG = FunctionSpace(mesh, "DG", 0)
V_vector = local_project(sqrt(inner(v, v)), DG).vector().get_local()
V_vector = local_project(sqrt(inner(v, v)), DG, local_rhs).vector().get_local()
h = mesh.hmin()

# Calculate flow rate at the inlet
Expand Down

0 comments on commit 6cfcc39

Please sign in to comment.