Skip to content

Linear

exponax.normalized.NormalizedLinearStepper ¤

Bases: BaseStepper

Source code in exponax/normalized/_linear.py
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
class NormalizedLinearStepper(BaseStepper):
    normalized_coefficients: tuple[float, ...]

    def __init__(
        self,
        num_spatial_dims: int,
        num_points: int,
        *,
        normalized_coefficients: tuple[float, ...] = (0.0, -0.5, 0.01),
    ):
        """
        Timestepper for d-dimensional (`d ∈ {1, 2, 3}`) linear PDEs on periodic
        boundary conditions with normalized dynamics.

        If the PDE in physical description is

            uₜ = ∑ᵢ aᵢ (∂ₓ)ⁱ u

        with `aᵢ` the coefficients on `domain_extent` `L` with time step `dt`,
        the normalized coefficients are

            αᵢ = (aᵢ Δt)/(Lⁱ)

        Important: note that the `domain_extent` is raised to the order of
        linear derivative `i`.

        One can also think of normalized dynamics, as a PDE on `domain_extent`
        `1.0` with time step `dt=1.0`.

        Take care of the signs!

        In the defaulf configuration of this timestepper, the PDE is an
        advection-diffusion equation with normalized advection of 0.5 and
        normalized diffusion of 0.01.

        **Arguments:**

        - `num_spatial_dims`: The number of spatial dimensions `d`.
        - `num_points`: The number of points `N` used to discretize the
            domain. This **includes** the left boundary point and **excludes**
            the right boundary point. In higher dimensions; the number of points
            in each dimension is the same. Hence, the total number of degrees of
            freedom is `Nᵈ`.
        - `normalized_coefficients`: The coefficients of the normalized
            dynamics. This must a tuple of floats. The length of the tuple
            defines the highest occuring linear derivative in the PDE.
        """
        self.normalized_coefficients = normalized_coefficients
        super().__init__(
            num_spatial_dims=num_spatial_dims,
            domain_extent=1.0,  # Derivative operator is just scaled with 2 * jnp.pi
            num_points=num_points,
            dt=1.0,
            num_channels=1,
            order=0,
        )

    def _build_linear_operator(self, derivative_operator: Array) -> Array:
        linear_operator = sum(
            jnp.sum(
                c * (derivative_operator) ** i,
                axis=0,
                keepdims=True,
            )
            for i, c in enumerate(self.normalized_coefficients)
        )
        return linear_operator

    def _build_nonlinear_fun(self, derivative_operator: Array):
        return ZeroNonlinearFun(
            num_spatial_dims=self.num_spatial_dims,
            num_points=self.num_points,
        )
__init__ ¤
__init__(
    num_spatial_dims: int,
    num_points: int,
    *,
    normalized_coefficients: tuple[float, ...] = (
        0.0,
        -0.5,
        0.01,
    )
)

Timestepper for d-dimensional (d ∈ {1, 2, 3}) linear PDEs on periodic boundary conditions with normalized dynamics.

If the PDE in physical description is

uₜ = ∑ᵢ aᵢ (∂ₓ)ⁱ u

with aᵢ the coefficients on domain_extent L with time step dt, the normalized coefficients are

αᵢ = (aᵢ Δt)/(Lⁱ)

Important: note that the domain_extent is raised to the order of linear derivative i.

One can also think of normalized dynamics, as a PDE on domain_extent 1.0 with time step dt=1.0.

Take care of the signs!

In the defaulf configuration of this timestepper, the PDE is an advection-diffusion equation with normalized advection of 0.5 and normalized diffusion of 0.01.

Arguments:

  • num_spatial_dims: The number of spatial dimensions d.
  • num_points: The number of points N used to discretize the domain. This includes the left boundary point and excludes the right boundary point. In higher dimensions; the number of points in each dimension is the same. Hence, the total number of degrees of freedom is Nᵈ.
  • normalized_coefficients: The coefficients of the normalized dynamics. This must a tuple of floats. The length of the tuple defines the highest occuring linear derivative in the PDE.
Source code in exponax/normalized/_linear.py
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
def __init__(
    self,
    num_spatial_dims: int,
    num_points: int,
    *,
    normalized_coefficients: tuple[float, ...] = (0.0, -0.5, 0.01),
):
    """
    Timestepper for d-dimensional (`d ∈ {1, 2, 3}`) linear PDEs on periodic
    boundary conditions with normalized dynamics.

    If the PDE in physical description is

        uₜ = ∑ᵢ aᵢ (∂ₓ)ⁱ u

    with `aᵢ` the coefficients on `domain_extent` `L` with time step `dt`,
    the normalized coefficients are

        αᵢ = (aᵢ Δt)/(Lⁱ)

    Important: note that the `domain_extent` is raised to the order of
    linear derivative `i`.

    One can also think of normalized dynamics, as a PDE on `domain_extent`
    `1.0` with time step `dt=1.0`.

    Take care of the signs!

    In the defaulf configuration of this timestepper, the PDE is an
    advection-diffusion equation with normalized advection of 0.5 and
    normalized diffusion of 0.01.

    **Arguments:**

    - `num_spatial_dims`: The number of spatial dimensions `d`.
    - `num_points`: The number of points `N` used to discretize the
        domain. This **includes** the left boundary point and **excludes**
        the right boundary point. In higher dimensions; the number of points
        in each dimension is the same. Hence, the total number of degrees of
        freedom is `Nᵈ`.
    - `normalized_coefficients`: The coefficients of the normalized
        dynamics. This must a tuple of floats. The length of the tuple
        defines the highest occuring linear derivative in the PDE.
    """
    self.normalized_coefficients = normalized_coefficients
    super().__init__(
        num_spatial_dims=num_spatial_dims,
        domain_extent=1.0,  # Derivative operator is just scaled with 2 * jnp.pi
        num_points=num_points,
        dt=1.0,
        num_channels=1,
        order=0,
    )
__call__ ¤
__call__(
    u: Float[Array, "C ... N"]
) -> Float[Array, "C ... N"]

Performs a check

Source code in exponax/_base_stepper.py
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
def __call__(
    self,
    u: Float[Array, "C ... N"],
) -> Float[Array, "C ... N"]:
    """
    Performs a check
    """
    expected_shape = (self.num_channels,) + spatial_shape(
        self.num_spatial_dims, self.num_points
    )
    if u.shape != expected_shape:
        raise ValueError(
            f"Expected shape {expected_shape}, got {u.shape}. For batched operation use `jax.vmap` on this function."
        )
    return self.step(u)