The PDE Definition Interface
While ODEs $u' = f(u,p,t)$ can be defined by a user-function f, for PDEs the function form can be different for every PDE. How many functions, and how many inputs? This can always change. The SciML ecosystem solves this problem by using ModelingToolkit.jl to define PDESystem, a high-level symbolic description of the PDE to be consumed by other packages.
The vision for the common PDE interface is that a user should only have to specify their PDE once, mathematically, and have instant access to everything as simple as a finite difference method with constant grid spacing, to something as complex as a distributed multi-GPU discontinuous Galerkin method.
The key to the common PDE interface is a separation of the symbolic handling from the numerical world. All of the discretizers should not "solve" the PDE, but instead be a conversion of the mathematical specification to a numerical problem. Preferably, the transformation should be to another ModelingToolkit.jl AbstractSystem via a symbolic_discretize dispatch, but in some cases this cannot be done or will not be performant. Thus in some cases, only a discretize definition is given to a AbstractSciMLProblem, with symbolic_discretize simply providing diagnostic or lower level information about the construction process.
These elementary problems, such as solving linear systems Ax=b, solving nonlinear systems f(x)=0, ODEs, etc. are all defined by SciMLBase.jl, which then numerical solvers can all target these common forms. Thus someone who works on linear solvers doesn't necessarily need to be working on a Discontinuous Galerkin or finite element library, but instead "linear solvers that are good for matrices A with properties ..." which are then accessible by every other discretization method in the common PDE interface.
Similar to the rest of the AbstractSystem types, transformation and analyses functions will allow for simplifying the PDE before solving it, and constructing block symbolic functions like Jacobians.
Constructors
ModelingToolkit.PDESystem — Typestruct PDESystem <: ModelingToolkit.AbstractSystemA system of partial differential equations.
Fields
- eqs: The equations which define the PDE.
- bcs: The boundary conditions.
- domain: The domain for the independent variables.
- ivs: The independent variables.
- dvs: The dependent variables.
- ps: The parameters.
- defaults: The default values to use when initial conditions and/or parameters are not supplied in- ODEProblem.
- connector_type: Type of the system.
- systems: The internal systems. These are required to have unique names.
- analytic: A vector of explicit symbolic expressions for the analytic solutions of each dependent variable. e.g.- analytic = [u(t, x) ~ a*sin(c*t) * cos(k*x)].
- analytic_func: A vector of functions for the analytic solutions of each dependent variable. Will be generated from- analyticif not provided. Should have the same argument signature as the variable, and a- psargument as the last argument, which takes an indexable of parameter values in the order you specified them in- ps. e.g.- analytic_func = [u(t, x) => (ps, t, x) -> ps[1]*sin(ps[2]*t) * cos(ps[3]*x)].
- name: The name of the system.
- description: A description of the system.
- metadata: Metadata for the system, to be used by downstream packages.
- gui_metadata: Metadata for MTK GUI.
Example
using ModelingToolkit
@parameters x t
@variables u(..)
Dxx = Differential(x)^2
Dtt = Differential(t)^2
Dt = Differential(t)
#2D PDE
C=1
eq  = Dtt(u(t,x)) ~ C^2*Dxx(u(t,x))
# Initial and boundary conditions
bcs = [u(t,0) ~ 0.,# for all t > 0
       u(t,1) ~ 0.,# for all t > 0
       u(0,x) ~ x*(1. - x), #for all 0 < x < 1
       Dt(u(0,x)) ~ 0. ] #for all  0 < x < 1]
# Space and time domains
domains = [t ∈ (0.0,1.0),
           x ∈ (0.0,1.0)]
@named pde_system = PDESystem(eq,bcs,domains,[t,x],[u])Domains (WIP)
Domains are specifying by saying indepvar in domain, where indepvar is a single or a collection of independent variables, and domain is the chosen domain type. A 2-tuple can be used to indicate an Interval. Thus forms for the indepvar can be like:
t ∈ (0.0, 1.0)
(t, x) ∈ UnitDisk()
[v, w, x, y, z] ∈ VectorUnitBall(5)Domain Types (WIP)
- Interval(a,b): Defines the domain of an interval from- ato- b(requires explicit import from- DomainSets.jl, but a 2-tuple can be used instead)
discretize and symbolic_discretize
The only functions which act on a PDESystem are the following:
- discretize(sys,discretizer): produces the outputted- AbstractSystemor- AbstractSciMLProblem.
- symbolic_discretize(sys,discretizer): produces a debugging symbolic description of the discretized problem.
Boundary Conditions (WIP)
Transformations
Analyses
Discretizer Ecosystem
NeuralPDE.jl: PhysicsInformedNN
NeuralPDE.jl defines the PhysicsInformedNN discretizer which uses a DiffEqFlux.jl neural network to solve the differential equation.
MethodOfLines.jl: MOLFiniteDifference (WIP)
MethodOfLines.jl defines the MOLFiniteDifference discretizer which performs a finite difference discretization using the DiffEqOperators.jl stencils. These stencils make use of NNLib.jl for fast operations on semi-linear domains.