Multigrid API
When using Decapodes to solve PDEs that use inverse Laplacians, you can use Multigrid to solve those linear systems.
API docs
CombinatorialSpaces.Multigrid.AbstractMultigridMode — Type
abstract type AbstractMultigridModeControls operator assembly in MultigridData construction.
DirectMode(): dualize and discretize at every level during the walk.GalerkinMode(): only dualize the finest mesh; derive coarse operators via the Galerkin condition Aₖ₊₁ = Rₖ Aₖ Pₖ.
CombinatorialSpaces.Multigrid.MeshTopology — Type
struct MeshTopologyBoundary maps of a simplicial 2-complex as plain arrays.
CombinatorialSpaces.Multigrid.MultigridData — Method
MultigridData(s, scheme, levels, op, steps; alg, mode)Construct a MultigridData from a base mesh using topology-only subdivision.
Keyword arguments
mode::AbstractMultigridMode = DirectMode(): controls operator assembly.alg = Circumcenter(): dualization algorithm.
Examples
md = MultigridData(s, BinarySubdivision(), 4, s -> ∇²(0, s), 3)
md = MultigridData(s, BinarySubdivision(), 4, s -> ∇²(0, s), 3; mode=GalerkinMode())CombinatorialSpaces.Multigrid._finalize_ops — Method
Direct mode: operators were collected coarsest-first; reverse to finest-first.
CombinatorialSpaces.Multigrid._finalize_ops — Method
Galerkin mode: dualize the finest mesh, derive coarse operators via Rₖ Aₖ Pₖ.
CombinatorialSpaces.Multigrid._init_ops — Method
Compute the coarsest-level operator before the walk begins.
CombinatorialSpaces.Multigrid._init_ops — Method
No operators needed before the walk in Galerkin mode.
CombinatorialSpaces.Multigrid._record_level_op! — Method
Dualize the current mesh and record its operator.
CombinatorialSpaces.Multigrid._record_level_op! — Method
No per-level operators in Galerkin mode.
CombinatorialSpaces.Multigrid.binary_subdivision_topo — Method
binary_subdivision_topo(topo::MeshTopology) -> (refined::MeshTopology, mat)Binary subdivision returning both topology and matrix.
CombinatorialSpaces.Multigrid.cubic_subdivision_topo — Method
cubic_subdivision_topo(topo::MeshTopology) -> (refined::MeshTopology, mat)Cubic subdivision returning both topology and matrix.
CombinatorialSpaces.Multigrid.full_multigrid — Function
Full multigrid: start at the coarsest grid and work up, applying μ-cycles at each level (μ=1 for V, μ=2 for W).
CombinatorialSpaces.Multigrid.make_restriction — Method
make_restriction(p::Transpose{<:Any, <:SparseMatrixCSC})Build restriction R = row_normalize(S) directly from the subdivision matrix S = parent(p), sharing its colptr and rowval arrays. Only allocates a new nzval vector and a temporary row-sum buffer.
CombinatorialSpaces.Multigrid.multigrid_vcycles — Function
Solve Ax=b with initial guess u for cycles V-cycles. alg is a Krylov.jl method (default cg).
CombinatorialSpaces.Multigrid.multigrid_wcycles — Function
W-cycle variant of multigrid_vcycles.
CombinatorialSpaces.Multigrid.propagate_points — Method
propagate_points(::BinarySubdivision, topo, coarse_points) -> fine_pointsInterpolate vertex positions for binary subdivision. Original vertices are copied; midpoints are averaged from edge endpoints.
CombinatorialSpaces.Multigrid.propagate_points — Method
propagate_points(::CubicSubdivision, topo, coarse_points) -> fine_pointsInterpolate vertex positions for cubic subdivision. Original vertices are copied; edge points at ⅓/⅔ positions and triangle centroids are computed from connectivity.
CombinatorialSpaces.Multigrid.propagate_points — Method
propagate_points(mat::SparseMatrixCSC, coarse_points) -> fine_pointsInterpolate via column-wise SpMV on the subdivision matrix. Exported for external callers that have a matrix but no scheme type; all internal paths use the scheme-dispatched methods above.
CombinatorialSpaces.Multigrid.refine — Method
refine(::BinarySubdivision, topo) -> MeshTopologyBinary (medial) topology refinement: split each edge at its midpoint, subdivide each triangle into 4. No subdivision matrix is constructed.
CombinatorialSpaces.Multigrid.refine — Method
refine(::CubicSubdivision, topo) -> MeshTopologyCubic topology refinement: two interior points per edge plus a centroid per triangle, subdividing each triangle into 9. No subdivision matrix is constructed.
CombinatorialSpaces.Multigrid.repeated_subdivision_maps — Method
repeated_subdivision_maps(k, ss, subdivider) -> Vector{PrimalGeometricMap}Apply k subdivisions via subdivider, returning geometric maps (meshes + matrices). Used by PrimalGeometricMapSeries.
CombinatorialSpaces.Multigrid.repeated_subdivisions — Method
repeated_subdivisions(k, ss, scheme) -> Vector{EmbeddedDeltaSet}Apply k subdivisions. Returns meshes only — no subdivision matrices are constructed.
CombinatorialSpaces.Multigrid.subdivision — Method
subdivision(s, scheme) -> EmbeddedDeltaSetSubdivide a mesh. Only refines topology and propagates points — no subdivision matrix is constructed.
CombinatorialSpaces.Multigrid.subdivision_map — Method
subdivision_map(s, scheme) -> PrimalGeometricMapSubdivide and return the geometric map (mesh + subdivision matrix).
CombinatorialSpaces.Multigrid.subdivision_matrix — Method
subdivision_matrix(::BinarySubdivision, topo) -> SparseMatrixCSCBuild the binary subdivision matrix (nvcoarse × nvfine) in direct CSC format. Identity block for original vertices, ½/½ averages for midpoints.
CombinatorialSpaces.Multigrid.subdivision_matrix — Method
subdivision_matrix(::CubicSubdivision, topo) -> SparseMatrixCSCBuild the cubic subdivision matrix (nvcoarse × nvfine) in direct CSC format. Identity block for original vertices, ⅓/⅔ weights for edge points, equal weights for centroids.
CombinatorialSpaces.Multigrid.topo_to_mesh — Method
topo_to_mesh(::Type{S}, topo::MeshTopology, points) -> SReconstitute an EmbeddedDeltaSet from a MeshTopology and point data.