scimba_torch.numerical_solvers.temporal_pde.time_discrete

Time-discrete numerical solvers for temporal PDEs.

Classes

ExplicitTimeDiscreteScheme(pde, projector[, ...])

An explicit time-discrete scheme for solving a differential equation.

TimeDiscreteAnagramProjector(pde, **kwargs)

A time-discrete natural gradient projector for solving temporal PDEs.

TimeDiscreteCollocationProjector(pde, **kwargs)

Implement the Galerkin-based nonlinear projection method.

TimeDiscreteImplicitNaturalGradientProjector(pde)

A time-discrete natural gradient for solving temporal PDEs.

TimeDiscreteNaturalGradientProjector(pde, ...)

A time-discrete natural gradient projector for solving temporal PDEs.

class TimeDiscreteCollocationProjector(pde, **kwargs)[source]

Bases: CollocationProjector

Implement the Galerkin-based nonlinear projection method.

This subclass implements the assembly method to assemble the input and output tensors for a specific nonlinear projection problem using the Galerkin approach. It computes the approximation of a nonlinear problem by sampling collocation points and evaluating the corresponding function values.

Parameters:
  • pde (TemporalPDE | KineticPDE) – The PDE model to solve.

  • **kwargs – Additional parameters for the projection, including collocation points and losses.

sample_all_vars(**kwargs)[source]

Samples all variables required for the projection.

Include the collocation and boundary points.

Parameters:

**kwargs (Any) – Additional keyword arguments.

Return type:

tuple[LabelTensor, ...]

Returns:

A tuple containing sampled collocation points and boundary data.

Raises:

ValueError – If the approximation space type is not recognized.

assembly_post_sampling(data, **kwargs)[source]

Assembles the I/Otensors for the nonlinear Galerkin projection problem.

This method samples collocation points and evaluates the corresponding function values from the approximation space and the right-hand side of the problem. It then returns the tensors representing the inputs and outputs of the projection.

Parameters:
  • data (tuple[LabelTensor, ...]) – tuple containing sampled collocation points and boundary data.

  • **kwargs – Additional keyword arguments.

Return type:

tuple[tuple[Tensor, ...], tuple[Tensor, ...]]

Returns:

A tuple of tensors representing the inputs (u) and outputs (f) of the projection problem.

Raises:

ValueError – If the approximation space type is not recognized.

class TimeDiscreteNaturalGradientProjector(pde, **kwargs)[source]

Bases: TimeDiscreteCollocationProjector

A time-discrete natural gradient projector for solving temporal PDEs.

Parameters:
  • pde (TemporalPDE | KineticPDE) – The PDE model to solve.

  • **kwargs – Additional parameters for the projection, including collocation points and losses.

class TimeDiscreteImplicitNaturalGradientProjector(pde, bc_type='strong', **kwargs)[source]

Bases: TimeDiscreteCollocationProjector

A time-discrete natural gradient for solving temporal PDEs.

Parameters:
  • pde (TemporalPDE | KineticPDE) – The PDE model to solve.

  • bc_type (str) – The way the boundary condition is handled; “strong” for strongly, “weak” for weakly

  • **kwargs – Additional parameters for the projection, including collocation points and losses.

class TimeDiscreteAnagramProjector(pde, **kwargs)[source]

Bases: TimeDiscreteCollocationProjector

A time-discrete natural gradient projector for solving temporal PDEs.

Parameters:
  • pde (TemporalPDE | KineticPDE) – The PDE model to solve.

  • **kwargs – Additional parameters for the projection, including collocation points and losses.

class ExplicitTimeDiscreteScheme(pde, projector, projector_init=None, **kwargs)[source]

Bases: object

An explicit time-discrete scheme for solving a differential equation.

Use linear and/or nonlinear spaces.

The class supports initialization of the model with a target function, computation of the right-hand side (RHS) from the model, and stepping through time using a projector.

Parameters:
initialization(**kwargs)[source]

Initializes the model by projecting the initial condition onto the model.

Parameters:

**kwargs (Any) – Additional parameters for the initialization, such as the number of epochs and collocation points.

abstract update(t, dt, **kwargs)[source]

Updates the model parameters using the time step and the chosen time scheme.

Parameters:
  • t (float) – The current time.

  • dt (float) – The time step.

  • **kwargs – Additional parameters for the update.

compute_relative_error_in_time(t, n_error=5000)[source]

Computes the relative error between the current and exact solution.

Parameters:
  • t (float) – The time at which the error is computed.

  • n_error (int) – The number of points used for computing the error. Default is 5_000.

Returns:

The L1, L2, and Linf errors.

Return type:

list

solve(dt=1e-05, final_time=0.1, sol_exact=None, **kwargs)[source]

Solves the full time-dependent problem, using time_step.

Parameters:
  • dt (float) – The time step.

  • final_time (float) – The final time.

  • sol_exact (Optional[Callable]) – The exact solution, if available.

  • **kwargs – Additional parameters for the time-stepping, such as the number of epochs, collocation points, and options for saving and plotting.