medypt.model
Base class for phase-field models.
Classes
Base class for phase-field models. |
Module Contents
- class medypt.model.ModelBase
Base class for phase-field models.
- opts: dict[str, Any]
A dictionary of numerical options. Initialized with default values and can be modified in subclasses. Contents include:
'has_thermal_noise'(bool): Whether to include thermal noise. Default isFalse.'rand_seed'(int): Random seed for thermal noise generation. Default is8347142.'quadr_deg'(int): Quadrature degree for numerical integration. Default is6.'petsc'(dict): A dictionary of PETSc solver options (see this for example). Default to use Newton nonlinear solver with MUMPS direct linear solver.'t_step_rtol'(float): Relative tolerance for time discretization error. Default is0.01.'dt_min_rescalar'(float): Minimum factor to reduce time step upon failure. Default is0.2.'dt_max_rescalar'(float): Maximum factor to increase time step upon success. Default is4.0.'dt_reducer'(float): Factor to reduce time step rescalar. Default is0.9.'max_successive_fail'(int): Maximum number of successive failures before stopping. Default is100.'min_dt'(float): Minimum time step size. Default is1e-9.'max_dt'(float): Maximum time step size. Default is10.0.'save_period'(int): Time step period for saving solution. Default is1.'log_file_name'(str): File name for logging evolution. Default isevolution.txt.'sol_file_name'(str): File name for saving solution. Default issolution.xdmf. If the suffix is not.xdmfor if no suffix, save solutions usingdolfinx.io.VTXWriterinto a folder with the given name.'verbose'(bool): Whether to print verbose output. Default isFalse.
- params: dict[str, Any]
A dictionary of physical parameters. Initialized to empty dictionary and should be set in subclasses.
- mesh_data: dolfinx.io.gmsh.MeshData
dolfinx.io.gmsh.MeshDataobject containing mesh and boundary tags. Initialized toNoneand should be set in subclasses.
- fields: dict[str, dolfinx.fem.Function]
A dictionary mapping field names to their
dolfinx.fem.Functionobjects. Initialized to empty dictionary and should be set in subclasses.
- load_mesh(comm: mpi4py.MPI.Comm, mesh: gmsh.model | str, **kwargs: Any)
Load the mesh data.
- Parameters:
comm (MPI.Comm) – MPI communicator.
mesh (gmsh.model | str) – Gmsh model or a file name with the
.mshor.xdmfformat.kwargs (dict[str, Any]) –
Additional keyword arguments passed to
create_mesh():mesh_dim (int): Geometric dimension of the mesh. Required when
meshis a gmsh model or a.mshfile.rank (int): Rank of the MPI process used for generating from gmsh model or reading from
.mshfiles. Required whenmeshis a gmsh model or a.mshfile.mesh_name (str): Name (identifier) of the mesh to read from the
.xdmffile. Required whenmeshis a.xdmffile.
- Returns:
None
- set_bcs(bcs: list[tuple[str, int | collections.abc.Callable[[Any], Any], dolfinx.fem.Constant | dolfinx.fem.Function | numpy.ndarray | collections.abc.Callable[[Any], Any]]])
Set essential and natural boundary conditions.
- Parameters:
bcs (list[tuple[str, int | Callable[[Any], Any], fem.Constant | fem.Function | np.ndarray | Callable[[Any], Any]]]) –
A list of tuples specifying boundary conditions. Each tuple is
(name, b, bc), where:name(str): The name of the field to which the boundary condition applies.b(int | Callable[[Any], Any]): The tag of or a callable defining the boundary. The tag is only for boundaries with codimension 1. The callable is only for Dirichlet boundary conditions and will have a calling signatureb(x), wherexis the spatial coordinate treated as a 2D numpy array with shape (geom_dim, num_points). It should return a boolean array of shape(num_points,)indicating whether each point is on the boundary.bc(Constant | Function | ndarray | Callable[[Any], Any]): A variable or callable defining the boundary condition. A variable is for Dirichlet boundary conditions, while a callable is for natural boundary conditions. The callable will have a calling signaturebc(fields), wherefieldsis treated as a tuple containing the splitted fields. It covers boundary conditions of the Neumann, Robin, and other general types implemented using the Nistche’s trick [1] [2]. The default is the zero-flux boundary condition.
Tip
One can make an evolving boundary condition by defining a global
dolfinx.fem.Constantordolfinx.fem.Functionobject and using it in the boundary-condition expression. One will then need to define an update function to update the value of the global object for a given time and pass that function tosolve().
- create_problem()
Create the finite element problem that is ready to be solved.
- solve(tspan: float, dt: float | None = None, ics: dict[str, collections.abc.Callable[[Any], Any] | dolfinx.fem.Function] | None = None, update: collections.abc.Callable[[float], None] = lambda t: ...) bool
Solve the time-dependent finite-element problem.
It uses adapative backward Euler time stepping algorithm for time integration. The time step is automatically controlled by the time discretization error specified in
opts. The advantage of backward Euler algorithm is its stability.- Parameters:
tspan (float) – Time span to solve over.
dt (float | None) – Initial time step. Defaults to
Noneto use the current time step.ics (dict[str, Callable[[Any], Any] | fem.Function] | None) – Initial conditions, defined as a dictionary mapping field names to callables that return the initial value for that field. The callables have a calling signature
ic(x), wherexis the spatial coordinate treated as a 2D numpy array with shape(geom_dim, num_points). Defaults toNoneto use the current values of the fields.update (Callable[[float], None]) – A callable to update any time-dependent parameters (must be defined a priori as global
dolfinx.fem.Constantordolfinx.fem.Function) at each time step. It takes the current time as the only argument. Defaults to a no-op function.
- Returns:
Trueif the solve completed successfully,Falseotherwise.- Return type:
bool