Phaseflow https://github.com/geo-fluid-dynamics/phaseflow-fenics

Documentation automatically generated by Sphinx

abstract_simulation.py provides an abstract class for unsteady adaptive FEM simulations.

We define a simulation as a sequence of time-dependent initial boundary value problems.

This module only provides an abstract base class for simulations. For an example of implementing a specific set of governing equations, see the abstract_phasechange_simulation module. For examples of implementing instantiable simulation classes, see e.g. cavity_melting_simulation and cavity_freezing_simulation.

This module essentially contains many lessons learned by a PhD student* applying FEniCS to time-dependent mixed FE problems with goal-oriented AMR. If you are not solving time-dependent problems, not using mixed finite elements, or not using goal-oriented AMR, then perhaps you do not need to consider much of this module.

The bulk of the difficult work is done by the fenics module; but fenics is mostly only a library for finite element spatial discretization of PDE’s. In that context, little attention is paid to the time discretization of initial value problem. Furthermore, we use some advanced features of fenics, namely mixed finite elements and goal-oriented adaptive mesh refinement (AMR), the latter of which is not among the better supported “core features” of fenics.

  • Alexander G. Zimmerman, AICES Graduate School, RWTH Aachen University
class phaseflow.abstract_simulation.AbstractSimulation(time_order=1, integration_measure=<MagicMock id='140089903650464'>, setup_solver=True)[source]

Bases: object

A class for time-dependent simulations with goal-oriented AMR using FeniCS

Attributes:
function_space
mesh
solution
time
timestep_size

Methods

adaptive_goal() Redefine this to return an adaptive goal.
advance() Move solutions backward in the queue to prepare for a new time step.
assign_initial_values() Set values of all solutions from self.initial_values().
boundary_conditions() Redefine this to return a list of fenics.DirichletBC.
coarse_mesh() Redefine this to return a fenics.Mesh.
deepcopy() Return an entire deep copy of self.
element() Redefine this to return a fenics.MixedElement.
governing_form() Redefine this to return a nonlinear variational form for the governing equations.
initial_mesh() Redefine this to refine the mesh before adaptive mesh refinement.
initial_values() Redefine this to return a fenics.Function containing the initial values.
load_newton_solution() When not using AMR, we can load a copy of the solution from the latest successful Newton iteration.
plot([solution_index, savefigs]) Plot the adaptive mesh and all parts of the mixed finite element solution.
print_constants() Print the names and values of all fenics.Constant attributes.
read_checkpoint(filepath) Read solutions and times from a checkpoint file.
reinit_solutions() Create the function space and solution functions for the current mesh and element.
reset_initial_guess() Set the values of the latest solution from the next latest solution.
save_newton_solution() When not using AMR, we can save a copy of the solution from the latest successful Newton iteration.
set_solution_on_subdomain(subdomain, values) Abuse fenics.DirichletBC to set values of a function on a subdomain.
setup_solver() Sometimes it is necessary to set up the solver again after breaking important references, e.g.
solve([goal_tolerance]) Solve the nonlinear variational problem.
time_discrete_terms() Apply first-order implicit Euler finite difference method.
write_checkpoint(filepath) Write solutions, times, and timestep sizes to a checkpoint file.
write_solution(file[, solution_index]) Write the solution to a file.
convert_checkpoints_to_xdmf_solution  
adaptive_goal()[source]

Redefine this to return an adaptive goal.

advance()[source]

Move solutions backward in the queue to prepare for a new time step. This is a separate method, since one may want to call solve multiple times before being satisfied with the solution state.

assign_initial_values()[source]

Set values of all solutions from self.initial_values().

boundary_conditions()[source]

Redefine this to return a list of fenics.DirichletBC.

coarse_mesh()[source]

Redefine this to return a fenics.Mesh.

deepcopy()[source]

Return an entire deep copy of self. For example, this is useful for checkpointing small problems in memory, or for running a batch of simulations with parameter changes.

element()[source]

Redefine this to return a fenics.MixedElement.

governing_form()[source]

Redefine this to return a nonlinear variational form for the governing equations.

initial_mesh()[source]

Redefine this to refine the mesh before adaptive mesh refinement.

initial_values()[source]

Redefine this to return a fenics.Function containing the initial values.

load_newton_solution()[source]

When not using AMR, we can load a copy of the solution from the latest successful Newton iteration. This can be useful, since a failed Newton iteration will blow up the solution, replacing it with garbage.

This will fail if the mesh has been changed by the adaptive solver and self.newton_solution has not been reinitialized with self.newton_solution = fenics.Function(self.function_space).

plot(solution_index=0, savefigs=False)[source]

Plot the adaptive mesh and all parts of the mixed finite element solution.

print_constants()[source]

Print the names and values of all fenics.Constant attributes. For example, this is useful for verifying that the correct parameters have been set.

read_checkpoint(filepath)[source]

Read solutions and times from a checkpoint file.

reinit_solutions()[source]

Create the function space and solution functions for the current mesh and element.

reset_initial_guess()[source]

Set the values of the latest solution from the next latest solution.

save_newton_solution()[source]

When not using AMR, we can save a copy of the solution from the latest successful Newton iteration. This can be useful, since a failed Newton iteration will blow up the solution, replacing it with garbage.

This will fail if the mesh has been changed by the adaptive solver and self.newton_solution has not been reinitialized with self.newton_solution = fenics.Function(self.function_space).

set_solution_on_subdomain(subdomain, values)[source]

Abuse fenics.DirichletBC to set values of a function on a subdomain.

Parameters:
subdomain

fenics.SubDomain

values

container of objects that would typically be passed to fenics.DirichletBC as the values of the boundary condition, one for each subspace of the mixed finite element solution space

setup_solver()[source]

Sometimes it is necessary to set up the solver again after breaking important references, e.g. after re-meshing.

solve(goal_tolerance=None)[source]

Solve the nonlinear variational problem. Optionally provide goal_tolerance to use the adaptive solver.

time_discrete_terms()[source]

Apply first-order implicit Euler finite difference method.

write_checkpoint(filepath)[source]

Write solutions, times, and timestep sizes to a checkpoint file.

write_solution(file, solution_index=0)[source]

Write the solution to a file. Parameters ———- file : fenics.XDMFFile

This method should have been called from within the context of the open file.
phaseflow.abstract_simulation.share_solver_parameters(share_to_parameters, share_from_parameters)[source]

FEniCS implements a setter for the solver parameters which does not allow us to

adaptive_solver.parameters[“nonlinear_variational_solver”] = solver.parameters

so we recursively catch the resulting KeyError exception and set all parameters.

abstract_phasechange_simulation.py provides an abstract class for phase-change simulations, including the natural and compositional convection of binary alloys.

class phaseflow.abstract_phasechange_simulation.AbstractPhaseChangeSimulation(time_order=1, integration_measure=<MagicMock id='140089903531176'>, setup_solver=True)[source]

Bases: phaseflow.abstract_simulation.AbstractSimulation

Implement the general model for phase-change coupled with natural and compositional convection.

This class is abstract, because an instantiable simulation still requires definitions for the mesh, initial values, and boundary conditions.

Attributes:
concentration_rayleigh_number
function_space
liquid_viscosity
liquidus_slope
mesh
prandtl_number
pressure_penalty_factor
pure_liquidus_temperature
regularization_central_temperature_offset
regularization_smoothing_parameter
schmidt_number
solid_viscosity
solution
stefan_number
temperature_rayleigh_number
time
timestep_size

Methods

adaptive_goal() Choose the solid area as the goal for AMR.
advance() Move solutions backward in the queue to prepare for a new time step.
assign_initial_values() Set values of all solutions from self.initial_values().
boundary_conditions() Redefine this to return a list of fenics.DirichletBC.
buoyancy(T, C) Extend the model from @cite{zimmerman2018monolithic} with a solute concentration.
coarse_mesh() Redefine this to return a fenics.Mesh.
coarsen([absolute_tolerances, …]) Re-mesh while preserving pointwise accuracy of solution variables.
deepcopy() Extends the parent deepcopy method with attributes for this derived class
element() Return a P1P2P1P1 element for the monolithic solution.
governing_form() Extend the model from @cite{zimmerman2018monolithic} with a solute concentration balance.
initial_mesh() Redefine this to refine the mesh before adaptive mesh refinement.
initial_values() Redefine this to return a fenics.Function containing the initial values.
load_newton_solution() When not using AMR, we can load a copy of the solution from the latest successful Newton iteration.
phi(T, C, T_m, m_L, delta_T, s) The regularized semi-phasefield.
plot([solution_index, savefigs]) Plot the adaptive mesh and all parts of the mixed finite element solution.
point_value_from_semi_phasefield(T, C) The semi-phasefield $phi$ sutiable for evaluation given a single $T$ and $C$.
print_constants() Print the names and values of all fenics.Constant attributes.
read_checkpoint(filepath) Read solutions and times from a checkpoint file.
reinit_solutions() Create the function space and solution functions for the current mesh and element.
reset_initial_guess() Set the values of the latest solution from the next latest solution.
save_newton_solution() When not using AMR, we can save a copy of the solution from the latest successful Newton iteration.
semi_phasefield(T, C) The semi-phasefield $phi$ given in UFL.
set_solution_on_subdomain(subdomain, values) Abuse fenics.DirichletBC to set values of a function on a subdomain.
setup_solver() Sometimes it is necessary to set up the solver again after breaking important references, e.g.
solve([goal_tolerance]) Solve the nonlinear variational problem.
solve_with_auto_regularization([…]) Catch solver failure and automatically over-regularize the problem, then successively return to desired regularization.
time_discrete_terms() Return the discrete time derivatives which are needed for the variational form.
write_checkpoint(filepath) Write solutions, times, and timestep sizes to a checkpoint file.
write_solution(file[, solution_index]) Write the solution to a file.
area_above_critical_phi_integrand  
convert_checkpoints_to_xdmf_solution  
solid_area_integrand  
solute_mass_integrand  
write_nonlinear_solver_table_header  
write_nonlinear_solver_table_row  
adaptive_goal()[source]

Choose the solid area as the goal for AMR.

buoyancy(T, C)[source]

Extend the model from @cite{zimmerman2018monolithic} with a solute concentration.

coarsen(absolute_tolerances=(0.01, 0.01, 0.01, 0.01, 0.01), maximum_refinement_cycles=6, circumradius_threshold=0.01)[source]

Re-mesh while preserving pointwise accuracy of solution variables.

deepcopy()[source]

Extends the parent deepcopy method with attributes for this derived class

element()[source]

Return a P1P2P1P1 element for the monolithic solution.

governing_form()[source]

Extend the model from @cite{zimmerman2018monolithic} with a solute concentration balance.

phi(T, C, T_m, m_L, delta_T, s)[source]

The regularized semi-phasefield.

point_value_from_semi_phasefield(T, C)[source]

The semi-phasefield $phi$ sutiable for evaluation given a single $T$ and $C$.

Maybe there is a way to evaluate the UFL expression rather than having to provide this redundant function.

semi_phasefield(T, C)[source]

The semi-phasefield $phi$ given in UFL.

solve_with_auto_regularization(goal_tolerance=None, max_regularization_threshold=4.0, max_attempts=16, enable_newton_solution_backup=False)[source]

Catch solver failure and automatically over-regularize the problem, then successively return to desired regularization.

If not using AMR, then the latest successful Newton solution can be saved/loaded to be more efficient with enable_newton_solution_backup = True.

time_discrete_terms()[source]

Return the discrete time derivatives which are needed for the variational form.

write_solution(file, solution_index=0)[source]

Write the solution to a file. Parameters ———- file : phaseflow.helpers.SolutionFile

This method should have been called from within the context of the open file.

abstract_heated_cavity_phasechange_simulation.py applies Phaseflow to benchmark problems on the unit square with hot and cold vertical walls, and adiabatic horizontal walls.

class phaseflow.abstract_heated_cavity_phasechange_simulation.AbstractHeatedCavityPhaseChangeSimulation(time_order=1, integration_measure=<MagicMock name='mock()' id='140089903778840'>, initial_uniform_gridsize=20, setup_solver=True)[source]

Bases: phaseflow.abstract_phasechange_simulation.AbstractPhaseChangeSimulation

Attributes:
concentration_rayleigh_number
function_space
liquid_viscosity
liquidus_slope
mesh
prandtl_number
pressure_penalty_factor
pure_liquidus_temperature
regularization_central_temperature_offset
regularization_smoothing_parameter
schmidt_number
solid_viscosity
solution
stefan_number
temperature_rayleigh_number
time
timestep_size

Methods

adaptive_goal() Choose the solid area as the goal for AMR.
advance() Move solutions backward in the queue to prepare for a new time step.
assign_initial_values() Set values of all solutions from self.initial_values().
boundary_conditions() Redefine this to return a list of fenics.DirichletBC.
buoyancy(T, C) Extend the model from @cite{zimmerman2018monolithic} with a solute concentration.
coarse_mesh() Redefine this to return a fenics.Mesh.
coarsen([absolute_tolerances, …]) Re-mesh while preserving pointwise accuracy of solution variables.
deepcopy() Extends the parent deepcopy method with attributes for this derived class
element() Return a P1P2P1P1 element for the monolithic solution.
governing_form() Extend the model from @cite{zimmerman2018monolithic} with a solute concentration balance.
initial_mesh() Redefine this to refine the mesh before adaptive mesh refinement.
initial_values() Redefine this to return a fenics.Function containing the initial values.
load_newton_solution() When not using AMR, we can load a copy of the solution from the latest successful Newton iteration.
phi(T, C, T_m, m_L, delta_T, s) The regularized semi-phasefield.
plot([solution_index, savefigs]) Plot the adaptive mesh and all parts of the mixed finite element solution.
point_value_from_semi_phasefield(T, C) The semi-phasefield $phi$ sutiable for evaluation given a single $T$ and $C$.
print_constants() Print the names and values of all fenics.Constant attributes.
read_checkpoint(filepath) Read solutions and times from a checkpoint file.
reinit_solutions() Create the function space and solution functions for the current mesh and element.
reset_initial_guess() Set the values of the latest solution from the next latest solution.
save_newton_solution() When not using AMR, we can save a copy of the solution from the latest successful Newton iteration.
semi_phasefield(T, C) The semi-phasefield $phi$ given in UFL.
set_solution_on_subdomain(subdomain, values) Abuse fenics.DirichletBC to set values of a function on a subdomain.
setup_solver() Sometimes it is necessary to set up the solver again after breaking important references, e.g.
solve([goal_tolerance]) Solve the nonlinear variational problem.
solve_with_auto_regularization([…]) Catch solver failure and automatically over-regularize the problem, then successively return to desired regularization.
time_discrete_terms() Return the discrete time derivatives which are needed for the variational form.
write_checkpoint(filepath) Write solutions, times, and timestep sizes to a checkpoint file.
write_solution(file[, solution_index]) Write the solution to a file.
area_above_critical_phi_integrand  
cold_wall_heat_flux_integrand  
convert_checkpoints_to_xdmf_solution  
solid_area_integrand  
solute_mass_integrand  
write_nonlinear_solver_table_header  
write_nonlinear_solver_table_row  
boundary_conditions()[source]

Redefine this to return a list of fenics.DirichletBC.

coarse_mesh()[source]

Redefine this to return a fenics.Mesh.

deepcopy()[source]

Extends the parent deepcopy method with attributes for this derived class

cavity_melting_simulation.py applies Phaseflow to the melting of a binary alloy in a square cavity.

class phaseflow.cavity_melting_simulation.CavityMeltingSimulation(initial_uniform_gridsize=1, time_order=1, integration_measure=<MagicMock name='mock()' id='140089903867776'>, setup_solver=True)[source]

Bases: phaseflow.abstract_heated_cavity_phasechange_simulation.AbstractHeatedCavityPhaseChangeSimulation

Attributes:
concentration_rayleigh_number
function_space
liquid_viscosity
liquidus_slope
mesh
prandtl_number
pressure_penalty_factor
pure_liquidus_temperature
regularization_central_temperature_offset
regularization_smoothing_parameter
schmidt_number
solid_viscosity
solution
stefan_number
temperature_rayleigh_number
time
timestep_size

Methods

adaptive_goal() Choose the solid area as the goal for AMR.
advance() Move solutions backward in the queue to prepare for a new time step.
assign_initial_values() Set values of all solutions from self.initial_values().
boundary_conditions() Redefine this to return a list of fenics.DirichletBC.
buoyancy(T, C) Extend the model from @cite{zimmerman2018monolithic} with a solute concentration.
coarse_mesh() Redefine this to return a fenics.Mesh.
coarsen([absolute_tolerances, …]) Re-mesh while preserving pointwise accuracy of solution variables.
deepcopy() Extends the parent deepcopy method with attributes for this derived class
element() Return a P1P2P1P1 element for the monolithic solution.
governing_form() Extend the model from @cite{zimmerman2018monolithic} with a solute concentration balance.
initial_mesh() Redefine this to refine the mesh before adaptive mesh refinement.
initial_values() Redefine this to return a fenics.Function containing the initial values.
load_newton_solution() When not using AMR, we can load a copy of the solution from the latest successful Newton iteration.
phi(T, C, T_m, m_L, delta_T, s) The regularized semi-phasefield.
plot([solution_index, savefigs]) Plot the adaptive mesh and all parts of the mixed finite element solution.
point_value_from_semi_phasefield(T, C) The semi-phasefield $phi$ sutiable for evaluation given a single $T$ and $C$.
print_constants() Print the names and values of all fenics.Constant attributes.
read_checkpoint(filepath) Read solutions and times from a checkpoint file.
reinit_solutions() Create the function space and solution functions for the current mesh and element.
reset_initial_guess() Set the values of the latest solution from the next latest solution.
save_newton_solution() When not using AMR, we can save a copy of the solution from the latest successful Newton iteration.
semi_phasefield(T, C) The semi-phasefield $phi$ given in UFL.
set_solution_on_subdomain(subdomain, values) Abuse fenics.DirichletBC to set values of a function on a subdomain.
setup_solver() Sometimes it is necessary to set up the solver again after breaking important references, e.g.
solve([goal_tolerance]) Solve the nonlinear variational problem.
solve_with_auto_regularization([…]) Catch solver failure and automatically over-regularize the problem, then successively return to desired regularization.
time_discrete_terms() Return the discrete time derivatives which are needed for the variational form.
write_checkpoint(filepath) Write solutions, times, and timestep sizes to a checkpoint file.
write_solution(file[, solution_index]) Write the solution to a file.
area_above_critical_phi_integrand  
cold_wall_heat_flux_integrand  
convert_checkpoints_to_xdmf_solution  
solid_area_integrand  
solute_mass_integrand  
write_nonlinear_solver_table_header  
write_nonlinear_solver_table_row  
deepcopy()[source]

Extends the parent deepcopy method with attributes for this derived class

initial_mesh()[source]

Redefine this to refine the mesh before adaptive mesh refinement.

initial_values()[source]

Redefine this to return a fenics.Function containing the initial values.

cavity_freezing_simulation.py applies Phaseflow to the freezing of a binary alloy in a square cavity.

class phaseflow.cavity_freezing_simulation.CavityFreezingSimulation(time_order=1, integration_measure=<MagicMock name='mock()' id='140089903855040'>, uniform_gridsize=20, setup_solver=True)[source]

Bases: phaseflow.abstract_heated_cavity_phasechange_simulation.AbstractHeatedCavityPhaseChangeSimulation

Attributes:
concentration_rayleigh_number
function_space
liquid_viscosity
liquidus_slope
mesh
prandtl_number
pressure_penalty_factor
pure_liquidus_temperature
regularization_central_temperature_offset
regularization_smoothing_parameter
schmidt_number
solid_viscosity
solution
stefan_number
temperature_rayleigh_number
time
timestep_size

Methods

adaptive_goal() Choose the solid area as the goal for AMR.
advance() Move solutions backward in the queue to prepare for a new time step.
assign_initial_values() Set values of all solutions from self.initial_values().
boundary_conditions() Redefine this to return a list of fenics.DirichletBC.
buoyancy(T, C) Extend the model from @cite{zimmerman2018monolithic} with a solute concentration.
coarse_mesh() Redefine this to return a fenics.Mesh.
coarsen([absolute_tolerances, …]) Re-mesh while preserving pointwise accuracy of solution variables.
deepcopy() Extends the parent deepcopy method with attributes for this derived class
element() Return a P1P2P1P1 element for the monolithic solution.
governing_form() Extend the model from @cite{zimmerman2018monolithic} with a solute concentration balance.
initial_mesh() Redefine this to refine the mesh before adaptive mesh refinement.
initial_values() Redefine this to return a fenics.Function containing the initial values.
load_newton_solution() When not using AMR, we can load a copy of the solution from the latest successful Newton iteration.
phi(T, C, T_m, m_L, delta_T, s) The regularized semi-phasefield.
plot([solution_index, savefigs]) Plot the adaptive mesh and all parts of the mixed finite element solution.
point_value_from_semi_phasefield(T, C) The semi-phasefield $phi$ sutiable for evaluation given a single $T$ and $C$.
print_constants() Print the names and values of all fenics.Constant attributes.
read_checkpoint(filepath) Read solutions and times from a checkpoint file.
reinit_solutions() Create the function space and solution functions for the current mesh and element.
reset_initial_guess() Set the values of the latest solution from the next latest solution.
save_newton_solution() When not using AMR, we can save a copy of the solution from the latest successful Newton iteration.
semi_phasefield(T, C) The semi-phasefield $phi$ given in UFL.
set_solution_on_subdomain(subdomain, values) Abuse fenics.DirichletBC to set values of a function on a subdomain.
setup_solver() Sometimes it is necessary to set up the solver again after breaking important references, e.g.
solve([goal_tolerance]) Solve the nonlinear variational problem.
solve_with_auto_regularization([…]) Catch solver failure and automatically over-regularize the problem, then successively return to desired regularization.
time_discrete_terms() Return the discrete time derivatives which are needed for the variational form.
write_checkpoint(filepath) Write solutions, times, and timestep sizes to a checkpoint file.
write_solution(file[, solution_index]) Write the solution to a file.
area_above_critical_phi_integrand  
cold_wall_heat_flux_integrand  
convert_checkpoints_to_xdmf_solution  
run  
set_constant_concentration  
setup_freezing_problem  
solid_area_integrand  
solute_mass_integrand  
solve_steady_state_heat_driven_cavity  
write_nonlinear_solver_table_header  
write_nonlinear_solver_table_row  
write_results_table_header  
write_results_table_row  
initial_mesh()[source]

Redefine this to refine the mesh before adaptive mesh refinement.

initial_values()[source]

Redefine this to return a fenics.Function containing the initial values.

refinement.py contains routines related to local mesh refinement.

phaseflow.refinement.adapt_coarse_solution_to_fine_solution(scalar, coarse_solution, fine_solution, element, absolute_tolerance=0.01, maximum_refinement_cycles=6, circumradius_threshold=0.01)[source]

Refine the mesh of the coarse solution until the interpolation error tolerance is met.

plotting.py defines utility functions for plotting meshes and finite element functions.

phaseflow.plotting.plot(f)[source]

This patches fenics.plot which is incorrect for functions on refind 1D meshes.

See https://bitbucket.org/fenics-project/dolfin/issues/1029/plotting-1d-function-incorrectly-ignores

helpers.py contains a variety of patching code.

phaseflow.helpers.mkdir_p(pathstring)[source]

Make a directory if it doesn’t exist.

This is needed because open does not create directories.

Now this just calls the appropriate function from pathlib. Older versions were more complicated.

Parameters:
path : string

backward_difference_formulas.py implements some BDF’s for time discretization.

phaseflow.backward_difference_formulas.apply_backward_euler(Delta_t, u)[source]

Apply the backward Euler (fully implicit, first order) time discretization method.

phaseflow.backward_difference_formulas.apply_bdf2(Delta_t, u)[source]

Apply the Gear/BDF2 (fully implicit, second order) backward difference formula for time discretization.

Here we use the variable time step size formula given by equation (12) from

@article{eckert2004bdf2,
title={A BDF2 integration method with step size control for elasto-plasticity}, author={Eckert, S and Baaser, H and Gross, D and Scherf, O}, journal={Computational Mechanics}, volume={34}, number={5}, pages={377–386}, year={2004}, publisher={Springer}

}

which we interpreted in the context of finite difference time discretizations in a PDE, rather than as an ODE solver, by comparing to the constant time step size BDF2 formula in

@article{belhamadia2012enhanced,
title={An enhanced mathematical model for phase change problems with natural convection}, author={Belhamadia, YOUSSEF and Kane, ABDOULAYE S and Fortin, ANDR{‘E}}, journal={Int. J. Numer. Anal. Model}, volume={3}, number={2}, pages={192–206}, year={2012}

}

Perhaps we could find a reference which derives variable time step size BDF2 in the context of finite difference time discretization, rather than in the context of ODE solvers.

Indices and tables