Poisson equation with periodic boundary conditions
Solving the Poisson equation on a periodic domain.
Introduction
This example solves the Poisson equation on the square with periodic boundary conditions instead of the usual Dirichlet conditions. Periodicity means the solution on the left edge equals the solution on the right edge, and similarly for the top and bottom edges:
Periodic boundary conditions are fundamental in:
- Homogenization of composite materials (solving cell problems)
- Wave propagation in periodic media (photonic/phononic crystals)
- Computational physics (periodic unit cells)
Implementation
Mesh Creation
We build a fine triangular mesh on a grid of points, then scale it to the domain :
#include <Rodin/Solver.h> #include <Rodin/Geometry.h> #include <Rodin/Variational.h> using namespace Rodin; using namespace Rodin::Geometry; using namespace Rodin::Variational; const size_t n = 200; Mesh mesh; mesh = mesh.UniformGrid(Polytope::Type::Triangle, { n, n }); mesh.getConnectivity().compute(1, 2); mesh.scale(1.0 / (n - 1)); // Normalize to [0, 1]^2 mesh.scale(2 * M_PI); // Stretch to [0, 2π]^2
Source Term
We use a smooth source function:
RealFunction f = [](const Point& p) { return sin(p.x() + M_PI / 4) * cos(p.y() + M_PI / 4); };
Building the DOF Correspondence Map
The key step for periodic conditions is building an IndexMap that pairs the DOF indices on opposite boundaries. For a uniform grid, the DOF layout is row-major. We pair:
- Left ↔ Right: DOF with DOF , for each row .
- Bottom ↔ Top: DOF with DOF , for each column .
P1 vh(mesh); IndexMap<IndexSet> dofs; // Left ↔ Right (pair first DOF in each row with last DOF in same row) for (Index i = 0; i < vh.getSize(); i += n) dofs[i].insert(i + n - 1); // Bottom ↔ Top (pair DOF j in first row with DOF j in last row) for (Index i = 0; i < n; i++) dofs[i].insert(i + n * (n - 1));
Problem Assembly and Solving
The problem is assembled identically to a standard Poisson equation, except that PeriodicBC replaces DirichletBC:
TrialFunction u(vh); TestFunction v(vh); Problem poisson(u, v); poisson = Integral(Grad(u), Grad(v)) - Integral(f, v) + PeriodicBC(u, dofs); Solver::SparseLU(poisson).solve();
Saving Results
u.getSolution().save("Periodic.gf", IO::FileFormat::MFEM); mesh.save("Periodic.mesh", IO::FileFormat::MFEM);
Full Source Code
The complete source code is:
See Also
- Poisson equation — Standard Dirichlet version
- Boundary conditions guide
- Weak formulations