Skip to content

第21章 网格分区(Mesh Partitioning)与区域分解


21.1 学习目标

  • 理解网格分区(Mesh Partitioning)即区域分解(Domain Decomposition)的概念,及其在大规模并行有限元计算中的角色
  • 掌握 Gmsh 四种分区器的选择与配置:Metis(图分区,默认)、Simple(几何切片)、Shell(内置回退)、Scotch 和 Hilbert(空间填充曲线)
  • 学会通过 mesh::partition(N)plugin::run("SimplePartition") 执行分区,以及 ONELAB 参数实现交互式分区控制
  • 理解分区实体(Partition Entity)的概念,掌握 getPartitions()getParent() 查询分区标签与父子关系
  • 掌握分区网格的导出选项:拓扑构建、Ghost 单元、物理组自动创建、分离文件模式

21.2 核心概念说明

21.2.1 为什么需要网格分区

在并行有限元计算中,单个网格可能包含数亿个单元,必须拆分为多个子域(subdomain)分发到不同计算节点。网格分区(Mesh Partitioning)即区域分解(Domain Decomposition),其目标是:

  1. 负载均衡(Load Balancing):各分区包含近似相等的计算代价
  2. 最小通信(Minimal Communication):分区间边界面/边数尽可能少
  3. 拓扑完好(Topology Preservation):分区边界表示完好,便于施加边界条件
原始网格(单域)                     分区后(N 个子域)
┌─────────────────┐              ┌──────┬──────┬──────┐
│    整体网格      │  ======>    │Part 0│Part 1│Part 2│
└─────────────────┘  partition   └──────┴──────┴──────┘

21.2.2 分区实体(Partition Entity)

Gmsh 分区时不修改原始网格文件,而是创建新的分区实体(Partition Entity)——一种特殊的离散基本实体(discrete elementary entity)。每个分区实体额外存储:

  • 分区索引(Partition Index):属于哪个分区(从 0 开始)
  • 父实体标签(Parent Entity Tag):指向分区前的原始实体
原始模型(2 个正方形,分区前):
  Surface 1 (dim=2, tag=1)    Surface 2 (dim=2, tag=2)

3 分区后(启用拓扑):
  分区 0: 离散面 10 ---parent---> Surface 1 的一部分
  分区 1: 离散面 11 ---parent---> Surface 1 的一部分 + Surface 2 的一部分
  分区 2: 离散面 12 ---parent---> Surface 2 的一部分
         各分区的离散边 ---parent---> 原始边

这种设计使分区后的模型仍然是完整的 B-rep 模型(如果启用 Mesh.PartitionCreateTopology),可继续施加边界条件、定义物理组。

21.2.3 分区器对比

通过选项 Mesh.Partitioner 选择:

分区器 算法 适用场景
0 Simple 轴对齐几何切片 规则域、调试
1 Metis(默认) 图分区(多级递归二分/K-way) 质量最优,需 METIS 库
2 Shell 时延定向图分区 METIS 不可用时的内置回退
3 Scotch 图分区(PT-Scotch) 分布式分区
4 Hilbert 空间填充曲线(SFC) 超大规模,速度最快

Metis 细化控制选项

选项 默认 说明
Mesh.MetisAlgorithm 1 1=递归二分, 2=K-way
Mesh.MetisObjective 1 1=最小边割, 2=最小通信量
Mesh.PartitionTriWeight 1.0 三角形单元计算权重
Mesh.PartitionQuadWeight 1.0 四边形单元计算权重
Mesh.PartitionTetWeight 1.0 四面体单元计算权重
Mesh.PartitionHexWeight 1.0 六面体单元计算权重

通过调整各单元类型权重,可使分区器在负载均衡时考虑不同单元类型的计算代价差异(如六面体计算代价是四面体的 3 倍,则设 HexWeight=3.0)。

21.2.4 关键分区选项

选项 默认 说明
Mesh.PartitionCreateTopology 1 为分区实体构建完整 B-rep(面、边、点)
Mesh.PartitionCreateGhostCells 0 每个分区额外包含一层相邻分区的 Ghost 单元
Mesh.PartitionCreatePhysicals 0 自动为分区实体创建物理组
Mesh.PartitionOldStyleMsh2 0 MSH2 兼容(旧格式,不支持分区实体)
Mesh.PartitionSplitMeshFiles 0 每分区输出独立文件(t21_0.msh, t21_1.msh...)

Ghost 单元:用于重叠区域分解(如 Schwarz 方法)。每个分区在边界外延一层邻接单元,避免在边界处构造特殊通信模式。非重叠方法(如 PETSc)通常不需要。

SimplePartition 插件:通过 plugin::setNumber("SimplePartition", "NumSlicesX/Y/Z", N) 设置各方向切片数,总分区数 = Nx * Ny * Nz。基于坐标而非图算法,速度快但负载均衡质量较差。


21.3 C++ 代码逐段讲解

21.3.1 几何建模与网格生成

#include <set>
#include <cmath>
#include <iostream>
#include <gmsh.h>

int main(int argc, char **argv)
{
  gmsh::initialize(argc, argv);
  gmsh::model::add("t21");

  // 创建两个相邻正方形,通过 fragment 在公共边处缝合
  gmsh::model::occ::addRectangle(0, 0, 0, 1, 1, 1);
  gmsh::model::occ::addRectangle(1, 0, 0, 1, 1, 2);
  std::vector<std::pair<int, int>> ov;
  std::vector<std::vector<std::pair<int, int>>> ovv;
  gmsh::model::occ::fragment({{2, 1}}, {{2, 2}}, ov, ovv);
  gmsh::model::occ::synchronize();

addRectangle(x,y,z,dx,dy,tag) 创建二维矩形面。fragment 将两个正方形在相交处共形分割,确保后续网格在公共边上节点一一对应——这是分区前的重要准备。

  // 设置网格尺寸并生成网格
  std::vector<std::pair<int, int>> tmp;
  gmsh::model::getEntities(tmp, 0);
  gmsh::model::mesh::setSize(tmp, 0.05);

  gmsh::model::addPhysicalGroup(2, {1}, 100, "Left");
  gmsh::model::addPhysicalGroup(2, {2}, 200, "Right");
  gmsh::model::mesh::generate(2);

必须先有网格才能分区。getEntities(tmp, 0) 获取所有 0 维实体(点),mesh::setSize 统一设置尺寸 0.05,然后生成二维网格。

21.3.2 ONELAB 参数定义

ONELAB 是 Gmsh 的交互式参数框架,让用户在 GUI 中实时调整分区设置而无需重新编译:

  gmsh::onelab::set(R"( [
  {
    "type":"number",
    "name":"Parameters/0Mesh partitioner",
    "values":[0],
    "choices":[0, 1],
    "valueLabels":{"Metis":0, "SimplePartition":1}
  },
  {
    "type":"number",
    "name":"Parameters/1Number of partitions",
    "values":[3],
    "min":1, "max":256, "step":1
  },
  {
    "type":"number",
    "name":"Parameters/2Create partition topology (BRep)?",
    "values":[1], "choices":[0, 1]
  },
  {
    "type":"number",
    "name":"Parameters/3Create ghost cells?",
    "values":[0], "choices":[0, 1]
  },
  {
    "type":"number",
    "name":"Parameters/3Create new physical groups?",
    "values":[0], "choices":[0, 1]
  },
  {
    "type":"number",
    "name":"Parameters/3Write file to disk?",
    "values":[1], "choices":[0, 1]
  },
  {
    "type":"number",
    "name":"Parameters/4Write one file per partition?",
    "values":[0], "choices":[0, 1]
  }
  ] )");

7 个 ONELAB 参数覆盖了分区器选择、分区数量、拓扑构建、ghost 单元、物理组创建、文件写入和分离文件模式。前缀相同的参数(如 Parameters/3*)在 GUI 中自动分组为同一行。

21.3.3 分区执行 Lambda

  auto partitionMesh = []()
  {
    std::vector<double> n;

    // 读取 ONELAB 参数并设置 Gmsh 选项
    gmsh::onelab::getNumber("Parameters/1Number of partitions", n);
    int N = static_cast<int>(n[0]);

    gmsh::onelab::getNumber("Parameters/2Create partition topology (BRep)?", n);
    gmsh::option::setNumber("Mesh.PartitionCreateTopology", n[0]);

    gmsh::onelab::getNumber("Parameters/3Create ghost cells?", n);
    gmsh::option::setNumber("Mesh.PartitionCreateGhostCells", n[0]);

    gmsh::onelab::getNumber("Parameters/3Create new physical groups?", n);
    gmsh::option::setNumber("Mesh.PartitionCreatePhysicals", n[0]);

    gmsh::option::setNumber("Mesh.PartitionOldStyleMsh2", 0);

    gmsh::onelab::getNumber("Parameters/4Write one file per partition?", n);
    gmsh::option::setNumber("Mesh.PartitionSplitMeshFiles", n[0]);

partitionMesh 是一个 lambda 函数,封装整个分区流程。它从 ONELAB 读取用户参数,设置对应的 Gmsh 选项,然后执行分区。当用户修改 GUI 参数并触发 "check" 事件时,该函数被重新调用以重新分区。

    // 选择分区器并执行
    gmsh::onelab::getNumber("Parameters/0Mesh partitioner", n);
    if(n[0] == 0) {
      gmsh::model::mesh::partition(N);    // Metis 图分区
    } else {
      gmsh::plugin::setNumber("SimplePartition", "NumSlicesX", N);
      gmsh::plugin::setNumber("SimplePartition", "NumSlicesY", 1);
      gmsh::plugin::setNumber("SimplePartition", "NumSlicesZ", 1);
      gmsh::plugin::run("SimplePartition");   // 几何切片
    }

两种分区方式:mesh::partition(N) 调用 Metis 图算法(分析单元邻接图找最优割);plugin::run("SimplePartition") 沿 X 轴均匀切片(NumSlicesY/Z=1 使切片为一维)。

    // 写入磁盘
    gmsh::onelab::getNumber("Parameters/3Write file to disk?", n);
    if(n[0]) gmsh::write("t21.msh");

Mesh.PartitionSplitMeshFiles=1write("t21.msh") 会为每个分区生成独立文件:t21_0.msh, t21_1.msh 等。

21.3.4 分区实体信息查询

    // 遍历所有实体,打印分区信息
    std::vector<std::pair<int, int>> entities;
    gmsh::model::getEntities(entities);

    for(auto e : entities) {
      std::vector<int> partitions;
      gmsh::model::getPartitions(e.first, e.second, partitions);
      if(partitions.size()) {
        std::string type;
        gmsh::model::getType(e.first, e.second, type);
        std::cout << "Entity (" << e.first << "," << e.second << ") "
                  << "of type " << type << "\n";
        std::cout << " - Partition(s):";
        for(auto p : partitions) std::cout << " " << p;
        std::cout << "\n";

        int pdim, ptag;
        gmsh::model::getParent(e.first, e.second, pdim, ptag);
        std::cout << " - Parent: (" << pdim << "," << ptag << ")\n";

        std::vector<std::pair<int, int>> bnd;
        gmsh::model::getBoundary({e}, bnd);
        std::cout << " - Boundary:";
        for(auto b : bnd)
          std::cout << " (" << b.first << "," << b.second << ")";
        std::cout << "\n";
      }
    }
  };

三个关键 API: - getPartitions(dim, tag, out):返回实体所属分区索引(空 = 非分区实体) - getParent(dim, tag, parentDim, parentTag):追溯分区实体的父实体 - getBoundary({e}, bnd):获取实体边界(分区实体的边界也是分区实体)

输出示例(3 分区,启用拓扑):

Entity (2,10) of type Discrete Surface
 - Partition(s): 0
 - Parent: (2,1)
 - Boundary: (1,30) (1,31) (1,32) (1,33)

"Discrete Surface" 说明分区实体是离散面——分区操作创建的是离散实体,不存储参数化几何。

21.3.5 GUI 事件循环

  partitionMesh();  // 首次执行

  // ONELAB check 事件回调:参数变更后重新分区
  auto checkForEvent = [=]() -> bool {
    std::vector<std::string> action;
    gmsh::onelab::getString("ONELAB/Action", action);
    if(action.size() && action[0] == "check") {
      gmsh::onelab::setString("ONELAB/Action", {""});
      partitionMesh();
      gmsh::graphics::draw();
    }
    return true;
  };

  std::set<std::string> args(argv, argv + argc);
  if(!args.count("-nopopup")) {
    gmsh::fltk::initialize();
    while(gmsh::fltk::isAvailable() && checkForEvent())
      gmsh::fltk::wait();
  }

  gmsh::finalize();
  return 0;
}

ONELAB 交互机制:用户修改参数后,ONELAB/Action 被设为 "check"checkForEvent 检测到后调用 partitionMesh() 重新分区并刷新显示。支持 -nopopup 命令行参数跳过 GUI。


21.4 完整可运行代码

// ch21_mesh_partitioning.cpp — 网格分区(区域分解)用于并行有限元
// 编译: g++ -o ch21_mesh_partitioning ch21_mesh_partitioning.cpp -lgmsh
#include <set>
#include <cmath>
#include <iostream>
#include <gmsh.h>

int main(int argc, char **argv)
{
    gmsh::initialize(argc, argv);
    gmsh::model::add("ch21");

    // ---- 1. 创建几何:两个相邻正方形 ----
    gmsh::model::occ::addRectangle(0, 0, 0, 1, 1, 1);
    gmsh::model::occ::addRectangle(1, 0, 0, 1, 1, 2);

    // ---- 2. 共形布尔分片:缝合公共边 ----
    std::vector<std::pair<int, int>> ov;
    std::vector<std::vector<std::pair<int, int>>> ovv;
    gmsh::model::occ::fragment({{2, 1}}, {{2, 2}}, ov, ovv);
    gmsh::model::occ::synchronize();

    // ---- 3. 设置网格尺寸并生成 ----
    std::vector<std::pair<int, int>> tmp;
    gmsh::model::getEntities(tmp, 0);
    gmsh::model::mesh::setSize(tmp, 0.05);

    gmsh::model::addPhysicalGroup(2, {1}, 100, "Left");
    gmsh::model::addPhysicalGroup(2, {2}, 200, "Right");
    gmsh::model::mesh::generate(2);

    // ---- 4. ONELAB 参数定义 ----
    gmsh::onelab::set(R"( [
    {
      "type":"number",
      "name":"Parameters/0Mesh partitioner",
      "values":[0],
      "choices":[0, 1],
      "valueLabels":{"Metis":0, "SimplePartition":1}
    },
    {
      "type":"number",
      "name":"Parameters/1Number of partitions",
      "values":[3],
      "min":1, "max":256, "step":1
    },
    {
      "type":"number",
      "name":"Parameters/2Create partition topology (BRep)?",
      "values":[1], "choices":[0, 1]
    },
    {
      "type":"number",
      "name":"Parameters/3Create ghost cells?",
      "values":[0], "choices":[0, 1]
    },
    {
      "type":"number",
      "name":"Parameters/3Create new physical groups?",
      "values":[0], "choices":[0, 1]
    },
    {
      "type":"number",
      "name":"Parameters/3Write file to disk?",
      "values":[1], "choices":[0, 1]
    },
    {
      "type":"number",
      "name":"Parameters/4Write one file per partition?",
      "values":[0], "choices":[0, 1]
    }
    ] )");

    // ---- 5. 分区执行函数 ----
    auto partitionMesh = []()
    {
        std::vector<double> n;

        gmsh::onelab::getNumber("Parameters/1Number of partitions", n);
        int N = static_cast<int>(n[0]);

        gmsh::onelab::getNumber("Parameters/2Create partition topology (BRep)?", n);
        gmsh::option::setNumber("Mesh.PartitionCreateTopology", n[0]);

        gmsh::onelab::getNumber("Parameters/3Create ghost cells?", n);
        gmsh::option::setNumber("Mesh.PartitionCreateGhostCells", n[0]);

        gmsh::onelab::getNumber("Parameters/3Create new physical groups?", n);
        gmsh::option::setNumber("Mesh.PartitionCreatePhysicals", n[0]);

        gmsh::option::setNumber("Mesh.PartitionOldStyleMsh2", 0);

        gmsh::onelab::getNumber("Parameters/4Write one file per partition?", n);
        gmsh::option::setNumber("Mesh.PartitionSplitMeshFiles", n[0]);

        gmsh::onelab::getNumber("Parameters/0Mesh partitioner", n);
        if(n[0] == 0) {
            gmsh::model::mesh::partition(N);
        } else {
            gmsh::plugin::setNumber("SimplePartition", "NumSlicesX", N);
            gmsh::plugin::setNumber("SimplePartition", "NumSlicesY", 1);
            gmsh::plugin::setNumber("SimplePartition", "NumSlicesZ", 1);
            gmsh::plugin::run("SimplePartition");
        }

        gmsh::onelab::getNumber("Parameters/3Write file to disk?", n);
        if(n[0]) gmsh::write("t21.msh");

        // 查询并打印分区实体信息
        std::vector<std::pair<int, int>> entities;
        gmsh::model::getEntities(entities);
        for(auto e : entities) {
            std::vector<int> partitions;
            gmsh::model::getPartitions(e.first, e.second, partitions);
            if(partitions.size()) {
                std::string type;
                gmsh::model::getType(e.first, e.second, type);
                std::cout << "Entity (" << e.first << "," << e.second << ") "
                          << "of type " << type << "\n";
                std::cout << " - Partition(s):";
                for(auto p : partitions) std::cout << " " << p;
                std::cout << "\n";

                int pdim, ptag;
                gmsh::model::getParent(e.first, e.second, pdim, ptag);
                std::cout << " - Parent: (" << pdim << "," << ptag << ")\n";

                std::vector<std::pair<int, int>> bnd;
                gmsh::model::getBoundary({e}, bnd);
                std::cout << " - Boundary:";
                for(auto b : bnd)
                    std::cout << " (" << b.first << "," << b.second << ")";
                std::cout << "\n";
            }
        }
    };

    // ---- 6. 首次分区 ----
    partitionMesh();

    // ---- 7. GUI 事件循环 ----
    auto checkForEvent = [=]() -> bool {
        std::vector<std::string> action;
        gmsh::onelab::getString("ONELAB/Action", action);
        if(action.size() && action[0] == "check") {
            gmsh::onelab::setString("ONELAB/Action", {""});
            partitionMesh();
            gmsh::graphics::draw();
        }
        return true;
    };

    std::set<std::string> args(argv, argv + argc);
    if(!args.count("-nopopup")) {
        gmsh::fltk::initialize();
        while(gmsh::fltk::isAvailable() && checkForEvent())
            gmsh::fltk::wait();
    }

    gmsh::finalize();
    return 0;
}

21.5 关键 API 速查表

21.5.1 分区执行

API 说明
gmsh::model::mesh::partition(N) 使用当前分区器(Metis/Scotch/Shell/Hilbert)划分 N 个分区
gmsh::plugin::run("SimplePartition") 执行 SimplePartition 插件(需先设置 NumSlicesX/Y/Z)

21.5.2 分区实体查询

API 说明
gmsh::model::getPartitions(dim, tag, out) 获取实体所属分区索引列表(空 = 非分区实体)
gmsh::model::getParent(dim, tag, parentDim, parentTag) 获取分区实体的父实体维度和标签

21.5.3 分区选项

选项 默认 说明
Mesh.Partitioner 1 0=Simple, 1=Metis, 2=Shell, 3=Scotch, 4=Hilbert
Mesh.PartitionCreateTopology 1 为分区实体构建 B-rep
Mesh.PartitionCreateGhostCells 0 创建一层 Ghost 单元
Mesh.PartitionCreatePhysicals 0 自动为分区实体创建物理组
Mesh.PartitionOldStyleMsh2 0 MSH2 兼容模式
Mesh.PartitionSplitMeshFiles 0 每分区输出独立文件

21.5.4 Metis 精细控制

选项 默认 说明
Mesh.MetisAlgorithm 1 1=递归二分, 2=K-way
Mesh.MetisObjective 1 1=最小边割, 2=最小通信量
Mesh.PartitionTriWeight 1.0 三角形权重
Mesh.PartitionQuadWeight 1.0 四边形权重
Mesh.PartitionTetWeight 1.0 四面体权重
Mesh.PartitionHexWeight 1.0 六面体权重

21.5.5 SimplePartition 插件参数

API / 参数 说明
plugin::setNumber("SimplePartition", "NumSlicesX", N) X 方向切片数
plugin::setNumber("SimplePartition", "NumSlicesY", N) Y 方向切片数
plugin::setNumber("SimplePartition", "NumSlicesZ", N) Z 方向切片数

21.5.6 ONELAB 交互

API 说明
onelab::set(json) 批量定义 ONELAB 参数(类型、范围、默认值)
onelab::getNumber(name, out) 读取数值型参数
onelab::getString(name, out) 读取字符串型参数

21.6 注意事项与最佳实践

  1. 分区必须在网格生成后mesh::partition()plugin::run("SimplePartition") 都要求已有网格。正确顺序:几何建模 -> 网格生成 -> 分区。

  2. 分区后不重新生成网格:分区后调用 mesh::generate() 会清除分区信息。分区应是最后一步操作。

  3. 分区实体是离散实体:类型字符串以 "Discrete" 开头(如 "Discrete Surface")。它们直接存储网格数据,无参数化几何表示。

  4. SimplePartition 限制:仅支持轴对齐均匀切片,复杂几何优先使用 Metis。

  5. 负载均衡与单元权重:对于各向异性网格(如边界层),调整单元类型权重以真实反映计算代价——例如边界层棱柱单元计算量高,可适当提高 Mesh.PartitionPriWeight

  6. 与 MPI 求解器对接:Gmsh 分区格式与大多数 MPI 并行求解器兼容。工作流:Gmsh 生成网格并分区 -> 导出 MSH4 -> 求解器读取每个分区的节点/单元 -> 各 MPI 进程独立组装求解。

  7. 分区质量评估:可以通过 getPartitions() 收集各分区单元数,计算最大/最小单元数比值来评估负载均衡度。理想值趋近于 1.0。

  8. 扩展阅读:参考 .geo 版本 t21.geo;扩展教程 x1.cpp 包含更多分区实体的深入操作。