2.4.5 Custom Class for coupled PDEs

This example shows how to solve a set of coupled PDEs, the spatially coupled FitzHugh–Nagumo model, which is a simple model for the excitable dynamics of coupled Neurons:

\[\begin{split}\partial_t u &= \nabla^2 u + u (u - \alpha) (1 - u) + w \\ \partial_t w &= \epsilon u\end{split}\]

Here, \(\alpha\) denotes the external stimulus and \(\epsilon\) defines the recovery time scale. We implement this as a custom PDE class below.

pde coupled
  0%|          | 0/100.0 [00:00<?, ?it/s]
Initializing:   0%|          | 0/100.0 [00:00<?, ?it/s]
  0%|          | 0/100.0 [00:00<?, ?it/s]
  0%|          | 0.31/100.0 [00:00<00:43,  2.31it/s]
  1%|          | 1.16/100.0 [00:00<00:15,  6.58it/s]
  5%|▌         | 5.28/100.0 [00:00<00:06, 13.89it/s]
 14%|█▍        | 14.39/100.0 [00:00<00:04, 17.29it/s]
 28%|██▊       | 27.95/100.0 [00:01<00:03, 18.55it/s]
 44%|████▍     | 44.47/100.0 [00:02<00:02, 19.04it/s]
 63%|██████▎   | 62.6/100.0 [00:03<00:01, 19.34it/s]
 82%|████████▏ | 81.71/100.0 [00:04<00:00, 19.52it/s]
 82%|████████▏ | 81.71/100.0 [00:05<00:01, 16.09it/s]
100%|██████████| 100.0/100.0 [00:05<00:00, 19.69it/s]
100%|██████████| 100.0/100.0 [00:05<00:00, 19.69it/s]

from pde import FieldCollection, PDEBase, UnitGrid


class FitzhughNagumoPDE(PDEBase):
    """FitzHugh–Nagumo model with diffusive coupling."""

    def __init__(self, stimulus=0.5, τ=10, a=0, b=0, bc="auto_periodic_neumann"):
        super().__init__()
        self.bc = bc
        self.stimulus = stimulus
        self.τ = τ
        self.a = a
        self.b = b

    def evolution_rate(self, state, t=0):
        v, w = state  # membrane potential and recovery variable

        v_t = v.laplace(bc=self.bc) + v - v**3 / 3 - w + self.stimulus
        w_t = (v + self.a - self.b * w) / self.τ

        return FieldCollection([v_t, w_t])


grid = UnitGrid([32, 32])
state = FieldCollection.scalar_random_uniform(2, grid)

eq = FitzhughNagumoPDE()
result = eq.solve(state, t_range=100, dt=0.01)
result.plot()

Total running time of the script: (0 minutes 5.271 seconds)