Beispiel #1
0
r"""
Kuramoto-Sivashinsky - Using `PDE` class
========================================

This example implements a scalar PDE using the :class:`~pde.pdes.pde.PDE`. We here
consider the `Kuramoto–Sivashinsky equation
<https://en.wikipedia.org/wiki/Kuramoto–Sivashinsky_equation>`_, which for instance 
describes the dynamics of flame fronts:

.. math::
    \partial_t u = -\frac12 |\nabla u|^2 - \nabla^2 u - \nabla^4 u
"""

from pde import PDE, ScalarField, UnitGrid

grid = UnitGrid([32, 32])  # generate grid
state = ScalarField.random_uniform(grid)  # generate initial condition

eq = PDE({"u": "-gradient_squared(u) / 2 - laplace(u + laplace(u))"})  # define the pde
result = eq.solve(state, t_range=10, dt=0.01)
result.plot()
Beispiel #2
0
def test_pde_consts():
    """test using the consts argument in PDE"""
    field = ScalarField(grids.UnitGrid([3]), 1)

    eq = PDE({"a": "b"}, consts={"b": 2})
    np.testing.assert_allclose(eq.evolution_rate(field).data, 2)

    eq = PDE({"a": "b**2"}, consts={"b": field})
    np.testing.assert_allclose(eq.evolution_rate(field).data, field.data)

    eq = PDE({"a": "laplace(b)"}, consts={"b": field})
    np.testing.assert_allclose(eq.evolution_rate(field).data, 0)

    eq = PDE({"a": "laplace(b)"}, consts={"b": 3})
    with pytest.raises(Exception):
        eq.evolution_rate(field)

    eq = PDE({"a": "laplace(b)"}, consts={"b": field.data})
    with pytest.raises(Exception):
        eq.evolution_rate(field)
Beispiel #3
0
This example implements a complex PDE using the :class:`~pde.pdes.pde.PDE`. We here
chose the `Schrödinger equation <https://en.wikipedia.org/wiki/Schrödinger_equation>`_ 
without a spatial potential in non-dimensional form:

.. math::
    i \partial_t \psi = -\nabla^2 \psi
    
Note that the example imposes Neumann conditions at the wall, so the wave packet is
expected to reflect off the wall.
"""

from math import sqrt

from pde import PDE, CartesianGrid, MemoryStorage, ScalarField, plot_kymograph

grid = CartesianGrid([[0, 20]], 128, periodic=False)  # generate grid

# create a (normalized) wave packet with a certain form as an initial condition
initial_state = ScalarField.from_expression(
    grid, "exp(I * 5 * x) * exp(-(x - 10)**2)")
initial_state /= sqrt(initial_state.to_scalar("norm_squared").integral.real)

eq = PDE({"ψ": f"I * laplace(ψ)"})  # define the pde

# solve the pde and store intermediate data
storage = MemoryStorage()
eq.solve(initial_state, t_range=2.5, dt=1e-5, tracker=[storage.tracker(0.02)])

# visualize the results as a space-time plot
plot_kymograph(storage, scalar="norm_squared")
Beispiel #4
0
"""
Time-dependent boundary conditions
==================================

This example solves a simple diffusion equation in one dimensions with time-dependent
boundary conditions.
"""

from pde import PDE, CartesianGrid, MemoryStorage, ScalarField, plot_kymograph

grid = CartesianGrid([[0, 10]], [64])  # generate grid
state = ScalarField(grid)  # generate initial condition

eq = PDE({"c": "laplace(c)"}, bc={"value_expression": "sin(t)"})

storage = MemoryStorage()
eq.solve(state, t_range=20, dt=1e-4, tracker=storage.tracker(0.1))

# plot the trajectory as a space-time plot
plot_kymograph(storage)
Beispiel #5
0
def test_pde_critical_input():
    """test some wrong input and edge cases"""
    # test whether reserved symbols can be used as variables
    grid = grids.UnitGrid([4])
    eq = PDE({"E": 1})
    res = eq.solve(ScalarField(grid), t_range=2)
    np.testing.assert_allclose(res.data, 2)

    with pytest.raises(ValueError):
        PDE({"t": 1})

    eq = PDE({"u": 1})
    assert eq.expressions == {"u": "1.0"}
    with pytest.raises(ValueError):
        eq.evolution_rate(FieldCollection.scalar_random_uniform(2, grid))

    eq = PDE({"u": 1, "v": 2})
    assert eq.expressions == {"u": "1.0", "v": "2.0"}
    with pytest.raises(ValueError):
        eq.evolution_rate(ScalarField.random_uniform(grid))

    eq = PDE({"u": "a"})
    with pytest.raises(RuntimeError):
        eq.evolution_rate(ScalarField.random_uniform(grid))

    eq = PDE({"x": "x"})
    with pytest.raises(ValueError):
        eq.evolution_rate(ScalarField(grid))
Beispiel #6
0
    
Here, :math:`D_0` and :math:`D_1` are the respective diffusivity and the
parameters :math:`a` and :math:`b` are related to reaction rates.

Note that the same result can also be achieved with a 
:doc:`full implementation of a custom class <pde_brusselator_class>`, which
allows for more flexibility at the cost of code complexity.
"""

from pde import PDE, FieldCollection, PlotTracker, ScalarField, UnitGrid

# define the PDE
a, b = 1, 3
d0, d1 = 1, 0.1
eq = PDE(
    {
        "u": f"{d0} * laplace(u) + {a} - ({b} + 1) * u + u**2 * v",
        "v": f"{d1} * laplace(v) + {b} * u - u**2 * v",
    }
)

# initialize state
grid = UnitGrid([64, 64])
u = ScalarField(grid, a, label="Field $u$")
v = b / a + 0.1 * ScalarField.random_normal(grid, label="Field $v$")
state = FieldCollection([u, v])

# simulate the pde
tracker = PlotTracker(interval=1, plot_args={"vmin": 0, "vmax": 5})
sol = eq.solve(state, t_range=20, dt=1e-3, tracker=tracker)
Beispiel #7
0
def test_pde_wrong_input():
    """test some wrong input"""
    with pytest.raises(ValueError):
        PDE({"t": 1})
    with pytest.raises(ValueError):
        PDE({"E": 1})
    with pytest.raises(ValueError):
        PDE({"E": 1, "t": 0})

    grid = grids.UnitGrid([4])
    eq = PDE({"u": 1})
    assert eq.expressions == {"u": "1.0"}
    with pytest.raises(ValueError):
        eq.evolution_rate(FieldCollection.scalar_random_uniform(2, grid))

    eq = PDE({"u": 1, "v": 2})
    assert eq.expressions == {"u": "1.0", "v": "2.0"}
    with pytest.raises(ValueError):
        eq.evolution_rate(ScalarField.random_uniform(grid))

    eq = PDE({"u": "a"})
    with pytest.raises(RuntimeError):
        eq.evolution_rate(ScalarField.random_uniform(grid))
    (25+30*tt+273.15,And(20<tt,tt<=25)), \
    (175+273.15,And(25<tt,tt<=197.5)), \
    (175+4*tt+273.15,And(197.5<tt,tt<=202.5)), \
    (195+273.15,And(202.5<tt,tt<=233)), \
    (195+8*tt+273.15,And(233<tt,tt<=238)), \
    (235+273.15,And(238<tt,tt<=268.5)), \
    (235+4*tt+273.15,And(268.5<tt,tt<=273.5)), \
    (255+273.15,And(273.5<tt,tt<=344.5)), \
    (255-(230/91)*tt+273.15,And(344.5<tt,tt<=435.5)), \
    (25+273.15,And(435.5<tt,tt<=500))).subs(tt,t))'

grid = CartesianGrid([[0, x_max]], [int(xstep)])
state=ScalarField(grid=grid,data=T_origin+T_Kelvin)

eq = PDE(rhs={'c': 'a*laplace(c)'}, 
    bc={'value_expression': bc_T},
    consts={'a':alpha})

storage = MemoryStorage()
eq.solve(state, t_range=tmax, tracker=storage.tracker(tstep))

def kymograph_pic():
    plot_kymograph(storage)

def temp_curve():
    p=[]
    for i in range(int(tmax/tstep)):
        p.append(storage.data[i][-1])
    
    q=np.asarray(p)
    plt.plot(np.linspace(0,tmax,int(tmax/tstep),endpoint=1),q,label='Temperature')
Beispiel #9
0
    ((175+20/{v}*(tt-197.5/{v})+{T_Kelvin}),And(197.5/{v}<tt,tt<=202.5/{v})),\
    ((195+{T_Kelvin}),And(202.5/{v}<tt,tt<=233/{v})),\
    ((195+40/{v}*(tt-233/{v})+{T_Kelvin}),And(233/{v}<tt,tt<=238/{v})),\
    ((235+{T_Kelvin}),And(238/{v}<tt,tt<=268.5/{v})),\
    ((235+20/{v}*(tt-268.5/{v})+{T_Kelvin}),And(268.5/{v}<tt,tt<=273.5/{v})),\
    ((255+{T_Kelvin}),And(273.5/{v}<tt,tt<=344.5/{v})),\
    ((255-230/{v}*(tt-344.5/{v})+{T_Kelvin}),And(344.5/{v}<tt,tt<=435.5/{v})),\
    ((25+{T_Kelvin}),true)).subs(tt,t))'

grid = CartesianGrid([[0, x_max]], [int(xstep)])
state = ScalarField(grid=grid, data=T_origin + T_Kelvin)

eq = PDE(
    rhs={'c': f'{alpha}*laplace(c)'},
    bc={
        'type': 'mixed',
        'value': -1 * h / llambda,
        'const': bc_T
    },
)
#bc_ops={})

storage = MemoryStorage()
eq.solve(state, t_range=(tmin, tmax), dt=1e-3, tracker=storage.tracker(tstep))


def temp_curve():
    p = []
    for i in range(int((tmax - tmin) / tstep)):
        p.append(storage.data[i][-1])

    q = np.asarray(p)
Beispiel #10
0
==============================

This example implements a PDE that is only defined in one dimension.
Here, we chose the `Korteweg-de Vries equation
<https://en.wikipedia.org/wiki/Korteweg–de_Vries_equation>`_, given by

.. math::
    \partial_t \phi = 6 \phi \partial_x \phi - \partial_x^3 \phi
    
which we implement using the :class:`~pde.pdes.pde.PDE`.
"""

from math import pi

from pde import PDE, CartesianGrid, MemoryStorage, ScalarField, plot_kymograph

# initialize the equation and the space
eq = PDE(
    {"φ": "6 * φ * get_x(gradient(φ)) - laplace(get_x(gradient(φ)))"},
    user_funcs={"get_x": lambda arr: arr[0]},
)
grid = CartesianGrid([[0, 2 * pi]], [32], periodic=True)
state = ScalarField.from_expression(grid, "sin(x)")

# solve the equation and store the trajectory
storage = MemoryStorage()
eq.solve(state, t_range=3, tracker=storage.tracker(0.1))

# plot the trajectory as a space-time plot
plot_kymograph(storage)
"""
1D problem - Using `PDE` class
==============================

This example implements a PDE that is only defined in one dimension.
Here, we chose the `Korteweg-de Vries equation
<https://en.wikipedia.org/wiki/Korteweg–de_Vries_equation>`_, given by

.. math::
    \partial_t \phi = 6 \phi \partial_x \phi - \partial_x^3 \phi
    
which we implement using the :class:`~pde.pdes.pde.PDE`.
"""

from math import pi

from pde import PDE, CartesianGrid, MemoryStorage, ScalarField, plot_kymograph

# initialize the equation and the space
eq = PDE({"φ": "6 * φ * d_dx(φ) - laplace(d_dx(φ))"})
grid = CartesianGrid([[0, 2 * pi]], [32], periodic=True)
state = ScalarField.from_expression(grid, "sin(x)")

# solve the equation and store the trajectory
storage = MemoryStorage()
eq.solve(state, t_range=3, tracker=storage.tracker(0.1))

# plot the trajectory as a space-time plot
plot_kymograph(storage)
using the :class:`~pde.pdes.pde.PDE` class. In particular, we consider
:math:`D(x) = 1.01 + \tanh(x)`, which gives a low diffusivity on the left side of the
domain.

Note that the naive implementation,
:code:`PDE({"c": "divergence((1.01 + tanh(x)) * gradient(c))"})`, has numerical
instabilities. This is because two finite difference approximations are nested. To
arrive at a more stable numerical scheme, it is advisable to expand the divergence, 

.. math::
    \partial_t c = D \nabla^2 c + \nabla D . \nabla c
"""

from pde import PDE, CartesianGrid, MemoryStorage, ScalarField, plot_kymograph

# Expanded definition of the PDE
diffusivity = "1.01 + tanh(x)"
term_1 = f"({diffusivity}) * laplace(c)"
term_2 = f"dot(gradient({diffusivity}), gradient(c))"
eq = PDE({"c": f"{term_1} + {term_2}"}, bc={"value": 0})

grid = CartesianGrid([[-5, 5]], 64)  # generate grid
field = ScalarField(grid, 1)  # generate initial condition

storage = MemoryStorage()  # store intermediate information of the simulation
res = eq.solve(field, 100, dt=1e-3,
               tracker=storage.tracker(1))  # solve the PDE

plot_kymograph(storage)  # visualize the result in a space-time plot