2D Advection-Diffusion Model (FEniCSx / dolfinx 0.10)#

AdvectionDiffusion2D([configs])

An implementation of the two dimensional Advection Diffusion Model using FEniCSx (dolfinx 0.9+).

create_AdvectionDiffusion2D_model([...])

A simple interface to creating an instance of AdvectionDiffusion2D.

A module that provides implementation(s) of simulation models implemented in FEniCSx (dolfinx 0.9+).

generate_ad_mesh(domain_size=(1.0, 1.0), buildings=None, resolution=0.05, mesh_file=None)[source]#

Generate a 2-D advection-diffusion mesh via Gmsh.

The domain is a rectangle [0, Lx] x [0, Ly] with optional rectangular building obstacles cut out of the interior.

Parameters:
  • domain_size (tuple) – (Lx, Ly) dimensions of the rectangular domain. Default (1.0, 1.0).

  • buildings

    list of building dicts, each with keys:

    • "origin"(x0, y0) lower-left corner of the building.

    • "width" — width (x direction).

    • "height" — height (y direction).

    Pass None (default) for a plain rectangular domain.

  • resolution (float) – characteristic mesh element size. Default 0.05.

  • mesh_file – optional path (str or Path) to save the generated mesh as an XDMF file (<mesh_file>.xdmf + <mesh_file>.h5). If None the mesh is only kept in memory.

Returns:

a dolfinx.mesh.Mesh.

Raises:

ImportError – if gmsh or the dolfinx Gmsh I/O is not available.

Example:

mesh = generate_ad_mesh(
    domain_size=(1.0, 1.0),
    buildings=[
        {"origin": (0.25, 0.15), "width": 0.25, "height": 0.25},
        {"origin": (0.60, 0.60), "width": 0.15, "height": 0.25},
    ],
    resolution=0.05,
)
class AdvectionDiffusion2DConfigs(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', model_name='AD-DOLFINx-2D', screen_output_iter=1, file_output_iter=1, time_integration=<factory>, num_prognostic_variables=None, space_discretization=<factory>, mesh_configs=<factory>, gls_stab=True, dt=0.2, t_eps=1e-06, sources=<factory>)[source]#

Bases: TimeDependentModelConfigs

Configurations for the two dimensional Advection Diffusion Model.

Parameters:
  • verbose (bool) – a boolean flag to control verbosity of the object.

  • debug (bool) – a boolean flag that enables adding extra functionality in a debug mode

  • output_dir (str | Path) – the base directory where the output files will be saved.

  • model_name (str) – name of the model. Default is None.

  • screen_output_iter (int) – iteration interval for screen output. Default is 1. Note that this should be a positive integer to enforce proper effect.

  • file_out_iter – iteration interval for file output. Default is 1. Note that this should be a positive integer to enforce proper effect.

  • time_integration (dict | None) –

    dictionary holding time integration configurations:

    • ’scheme’: string specifying the time integration scheme (e.g., ‘RK4’, ‘RK45’, ‘BDF’, etc.). Default is None.

    • ’stepsize’: float specifying the time integration stepsize. Default is None.

    • ’adaptive’: bool specifying whether the time integration is adaptive. Default is False.

  • num_prognostic_variables (int | None) – number of prognostic variables in the model. Default is None. Must be a positive integer if not None.

  • space_discretization (dict | None) –

    dictionary holding space discretization configurations. Contains:

    • ’scheme’: string specifying the space discretization scheme (e.g., ‘FD’, ‘FE’, ‘BE’, etc.). Default is None.

  • mesh_configs (dict) –

    dictionary holding mesh configurations/settings. They include:

    • mesh_file: path to the mesh file. Default path is to included ad_20.xml file. For Gmsh-based generation this is used as the save path (optional).

    • n_mesh_refinement: number of mesh refinements. Default is 0.

    • mesh_elements: type of mesh elements. Default is Lagrange.

    • mesh_elements_degree: degree of mesh elements. Default is 2.

    • buildings: None (default) or list of {origin, width, height} dicts. When non-None, Gmsh is used to generate the mesh.

    • resolution: float mesh element size for Gmsh (default 0.05).

  • gls_stab (bool) – use Galerkin Least Squares stabilization. Default is True.

  • dt (float) – time step size. Default is 0.2.

  • t_eps (float) – time precision. Default is 1e-6.

  • sources (list) –

    list of source dicts, each describing one Gaussian contaminant source. Each dict must contain:

    • 'center': (cx, cy) — 2-D source centre coordinates (required)

    • 'amplitude': float cap applied to the Gaussian (default 0.5)

    • 'width': Gaussian decay coefficient (default 100.0)

    Defaults to the canonical two-source layout used in published results.

mesh_configs: dict#
gls_stab: bool#
dt: float#
t_eps: float#
sources: list#
__init__(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', model_name='AD-DOLFINx-2D', screen_output_iter=1, file_output_iter=1, time_integration=<factory>, num_prognostic_variables=None, space_discretization=<factory>, mesh_configs=<factory>, gls_stab=True, dt=0.2, t_eps=1e-06, sources=<factory>)#
class AdvectionDiffusion2D(configs=None)[source]#

Bases: TimeDependentModel

An implementation of the two dimensional Advection Diffusion Model using FEniCSx (dolfinx 0.9+).

Parameters:

configs (AdvectionDiffusion2DConfigs | dict | None) – object holding model configurations/settings.

__init__(configs=None)[source]#
validate_configurations(configs, raise_for_invalid=True)[source]#

Check the passed configurations and make sure they are conformable with each other, and with current configurations once combined.

Note

Here only the locally-defined configurations in AdvectionDiffusion2DConfigs are validated. Finally, super classes validators are called.

Parameters:
Returns:

flag indicating whether passed configurations dictionary is valid or not

Raises:
Return type:

bool

parameter_vector(init_val=0, return_np=True)[source]#

Create an instance of model parameter vector.

Parameters:
  • init_val (float) – value assigned to entries of the parameter vector

  • return_np (bool) – return numpy array (True) or PETSc Vec (False)

Returns:

1-D numpy array or PETSc.Vec

state_vector(init_val=0, return_np=True)[source]#

Create an instance of model state vector.

Parameters:
  • init_val (float) – value assigned to entries of the state vector

  • return_np (bool) – return numpy array (True) or PETSc Vec (False)

Returns:

1-D numpy array or PETSc.Vec

is_state_vector(state)[source]#

Test whether the passed state vector is valid or not

is_parameter_vector(parameter)[source]#

Test whether the passed parameter vector is valid or not

integrate_state(state, tspan, checkpoints=None, return_np=True, verbose=False)[source]#

Simulate/integrate the model starting from the initial state over the passed checkpoints.

Parameters:
  • state – data structure holding the initial model state

  • tspan(t0, tf) iterable specifying the integration window

  • checkpoints – times at which to store the computed solution. If None (default), use [t0, t0+dt, ..., tf].

  • return_np (bool) – return numpy arrays (True, default) or PETSc Vecs

  • verbose (bool) – output progress to screen

Returns:

(checkpoints, trajectory) — a sorted array of times and a list of state snapshots at those times.

solve_forward_sub(ic, t_init, t_final, return_np=True)[source]#

Propagate the initial condition ic over the interval [t_init, t_final].

Parameters:
  • ic – initial condition (numpy array or PETSc.Vec)

  • t_init (float) – start time

  • t_final (float) – end time

  • return_np (bool) – return numpy array (True) or PETSc.Vec (False)

Returns:

state at t_final

Jacobian_T_matvec(state, eval_at_t=None, eval_at=None, dt=None, return_np=True)[source]#

Evaluate and return the product of the Jacobian transposed (adjoint operator) by a model state. This model is linear so the TLM is the model itself.

Parameters:
  • state – state to multiply the Jacobian transpose by

  • eval_at_t – time at which the Jacobian is evaluated (ignored, model is linear)

  • eval_at – state around which the Jacobian is evaluated (ignored)

  • dt (float) – the step size

  • return_np (bool) – return numpy array or PETSc.Vec

Returns:

product J^T * state

solve_adjoint_sub(adjoint, t_final, t_init, return_np=True)[source]#

Propagate the adjoint/sensitivity vector backward over [t_init, t_final].

Parameters:
  • adjoint – adjoint vector (numpy array or PETSc.Vec)

  • t_final (float) – end time

  • t_init (float) – start time

  • return_np (bool) – return numpy array (True) or PETSc.Vec (False)

Returns:

adjoint state at t_init

create_initial_condition(sources=None, return_np=True)[source]#

Create the initial condition vector.

Each source is a Gaussian bump capped at amplitude:

s_i(x) = min(amplitude_i, exp(-width_i * ((x-cx_i)^2 + (y-cy_i)^2)))

Multiple sources are combined with equal weights:

ic(x) = (1/n) * sum_i s_i(x)
Parameters:
  • sources

    None to use the model’s default source layout (controlled by configurations.num_sources), or a non-empty list of dicts, one per source. Each dict must contain:

    • 'center': (cx, cy) — 2-D source centre coordinates (required)

    • 'amplitude': float cap applied before summing (default 0.5)

    • 'width': Gaussian decay coefficient (default 100.0)

  • return_np (bool) – return numpy array (True) or PETSc.Vec (False)

Returns:

initial condition as numpy array or PETSc.Vec

Raises:
  • TypeError – if sources is not None or a non-empty list of dicts

  • ValueError – if any source dict is missing the required 'center' key, or if configurations.num_sources is not 1 or 2 when sources is None

get_model_grid()[source]#

Return a numpy array of the DOF coordinates, shape (ndofs, 2).

plot_domain(ax=None, clean_xy_ticks=True, xlim=None, ylim=None, title=None, show_mesh=True, mesh_linewidth=1.0, mesh_linecolor='#000C67', mesh_alpha=0.5, mesh_zorder=10, triangles_lw=0.01, boxes_coordinates=None, boxes_texts=None, boxes_hatch=None, boxes_zorder=5, sensors_coordinates=None, sensors_radius=15, annotate_sensors=False, sensors_alpha=0.75, sensors_zorder=10, prediction_coordinates=None, prediction_radius=2, prediction_alpha=0.85, prediction_zorder=15, fontsize=18, show_axis_frame=False, save_to_file=None, verbose=False)[source]#

Plot the domain, with optional buildings, sensor locations, and prediction coordinates.

Parameters:

ax – matplotlib axis or None (a new figure is created).

Returns:

matplotlib axis

plot_state(vec, ax=None, clean_xy_ticks=True, xlim=None, ylim=None, cmap='jet', alpha=0.75, add_colorbar=True, triangles_lw=0.01, show_mesh=True, mesh_linewidth=1.0, mesh_linecolor='#000C67', mesh_alpha=0.25, mesh_zorder=10, boxes_coordinates=None, boxes_texts=None, boxes_hatch=None, boxes_zorder=5, sensors_coordinates=None, sensors_radius=15, annotate_sensors=False, sensors_alpha=0.75, sensors_zorder=10, prediction_coordinates=None, prediction_radius=2, prediction_alpha=0.85, prediction_zorder=15, title=None, fontsize=18, show_axis_frame=True, save_to_file=None, verbose=False)[source]#

2D tricontourf plot of the passed state vector.

Parameters:

vec – model state (numpy array or PETSc.Vec)

Returns:

matplotlib axis

plot_velocity_field(ax=None, show_mesh=True, mesh_alpha=0.1, show_axis_frame=False, boxes_coordinates=None, boxes_texts=None, boxes_hatch=None, sensors_coordinates=None, sensors_radius=15, annotate_sensors=False, fontsize=18, title=None, save_to_file=None, verbose=False)[source]#

Plot the velocity field using a quiver plot.

The velocity DOF values and coordinates are extracted directly from the dolfinx Function — no temporary dolfin figure is needed.

Returns:

matplotlib axis

animate_state(state, tspan, checkpoints=None, clean_xy_ticks=True, xlim=None, ylim=None, cmap='jet', alpha=0.75, add_colorbar=True, triangles_lw=0.01, show_mesh=True, mesh_linewidth=1.0, mesh_linecolor='#000C67', mesh_alpha=0.25, mesh_zorder=10, boxes_hatch=None, boxes_zorder=5, sensors_coordinates=None, sensors_radius=15, annotate_sensors=False, sensors_alpha=0.75, sensors_zorder=10, prediction_coordinates=None, prediction_radius=2, prediction_alpha=0.85, prediction_zorder=15, fontsize=18, show_axis_frame=True, show_times=True, framerate=3, dpi=150, output_dir=None, video_file_name='AD_Simulation.mp4', return_axes=False)[source]#

Loop over plot_state() to create an animated video of the simulated trajectory.

Frames are piped directly to ffmpeg via matplotlib.animation.FFMpegWriter — no intermediate PNG files are written to disk.

Note

Requires ffmpeg to be installed and on the PATH.

Parameters:
  • show_times – set the simulation time as the frame title.

  • framerate – frames per second (positive integer).

  • dpi – resolution used when rendering each frame (default 150).

  • output_dir – directory in which to save the video. Defaults to the model’s output_dir config, then os.getcwd().

  • video_file_name – filename for the .mp4 output.

  • return_axes – if True return a {time: ax} dict of per-frame axes (each on its own persistent figure) instead of None.

Returns:

None, or dict of axes when return_axes=True.

Raises:

RuntimeError – if ffmpeg is not available or exits with an error.

property true_initial_condition#

The default (ground-truth) initial condition.

property state_dof#

Function space corresponding to the model state.

property parameter_dof#

Function space corresponding to the model parameter.

property mesh_configs#

Retrieve the model mesh configurations.

property gls_stab#

Retrieve the GLS stabilization flag.

property boxes_coordinates#

Automatically identify boxes coordinates based on the mesh file name.

property boxes_texts#

Default list (or None) of names associated with the boxes.

property dt#

Retrieve the model time step size.

property t_eps#

Retrieve the model time epsilon tolerance.

property t_precision#

Evaluate time precision for rounding to be used with np.round().

create_AdvectionDiffusion2D_model(mesh_file='/Users/attia/AHMED_HOME/Research/Projects/Software/PyOED/Branches/main/pyoed/doc/source/../../pyoed/models/simulation_models/fenics_models/meshes/ad_20.xml', mesh_elements='Lagrange', mesh_elements_degree=2, n_mesh_refinement=0, gls_stab=True, dt=0.2, output_dir='./_PYOED_RESULTS_')[source]#

A simple interface to creating an instance of AdvectionDiffusion2D.

Parameters:
  • mesh_file – path to the mesh file (XML or XDMF).

  • mesh_elements (str) – finite element family name.

  • mesh_elements_degree (int) – polynomial degree.

  • n_mesh_refinement (int) – number of uniform mesh refinements.

  • gls_stab (bool) – use GLS stabilization.

  • dt (float) – time step.

  • output_dir – output directory for plots/results.

Returns:

an AdvectionDiffusion2D instance.