AMCAX Kernel
Geometry kernel for CAD/CAE/CAM
九韶内核 1.0.0.0
载入中...
搜索中...
未找到
网格剖分示例

概述

本教程提供了使用网格剖分的基本用法,导入 STEP 格式或者 STL 格式文件,通过 json 文件传入网格剖分控制参数,自动生成体网格。

输入形状与输出网格预览

左图为导入的CAD模型,右图为最终生成的网格模型。

左图为导入的STL格式模型,右图为最终生成的网格模型。

命名空间

using namespace AMCAX::NextMesh;
AMCAX NextMesh 模块提供的所有接口所在的命名空间。
定义 misc.docu:42

准备工作

设置log格式,提供cad文件和网格生成控制json文件。

const std::string cadstr = "./data/qiuwotaojian.stp";
const std::string jsonPath = "./data/mesh_para-mesh000.json";
static AMCAX_API void InitLogger(const std::string &logFileName="", const std::string &logPattern="[%Y-%m-%dT%X%z] [%l] %v", const int64_t logFileSize=20 *1024 *1024, const int maxLogFiles=100)
设置日志记录器(默认情况下日志功能是关闭的,除非调用该函数)

Json示例

通过设置以下剖分参数控制网格生成:待剖分实体维度MeshingDim;网格单元最大最小尺寸MaxSize/MinSize(绝对尺寸),可用几何模型boundingbox对角线长度bblen乘以一定比例r,比如[0.01, 0.1]*bblen;网格增长率GrowthRate;待剖分实体Tag序列SelectedEntities;曲率因子CurvatureFactor;网格剖分类型ElementType;网格剖分方法MeshingMethod;并行线程数ThreadNum。

{
"MeshSize": [
{
"CurvatureFactor": 0.6,
"GrowthRate": 1.0,
"MaxSize": 0.5,
"MeshingDim": 3,
"MinSize": 0.1,
"ElementType": "Tri",
"MeshingMethod": "AFTDelaunay",
"SelectedEntities": []
}
],
"ThreadNum": 1
}

注意:当输入 STL 文件时,目前 “MeshingMethod” 只能使用 AFT。另外,当前尺寸设置对输入 STL 模型无效。

是否必需 取值范围 默认值 备注
MeshSize CurvatureFactore [0, 1.0] 0 0表示关闭;1.0表示整圆剖分成3段
MeshSize GrowthRate [1.0, inf] 1.0
MeshSize MaxSize [eps, inf] inf 绝对尺寸
MeshSize MinSize [0, inf] 0 绝对尺寸,MinSize<=MaxSize
MeshSize MeshingDim [1,2,3] ——
MeshSize SelectedEntities MeshingDim维度下实体tag的子集 —— 取空值表示指定维度下的所有实体
MeshSize ElementType [Tri, Quad, Mixed, RTri] Tri 只在MeshingDim=2时有效,其余忽略;Tri 三角形;Quad 全四边形;RTri 直角三角形;Mixed 三角形四边形混合网格
MeshSize MeshingMethod [Transfinite,Delaunay,AFTDelaunay,AFT] ElementType="Tri"时默认"AFTDelaunay" 只在MeshingDim=2时有效,其余忽略;Transfinite方法生成结构化网格,仅对含2/3/4边的face有效;Delaunay 德劳内方法生成三角形;AFTDelaunay 波前与德劳内结合方法生成三角形

Json示例-平面边界层

以二维翼型生成平面边界层混合网格为例。

{
"BoundaryLayer": [
{
"DistributionType": "Separation",
"EntityType": "Edge",
"EFPairs": [
{
"E": 5,
"F": 1
},
{
"E": 6,
"F": 1
}
],
"Basic": {
"GrowthRate": {
"Constant": 1.2
},
"InitLayerHeight": {
"Absolute": 0.01
},
"NLayers": 20,
"ProblematicAreasTreatment": "STOP"
}
}
],
"MeshSize": [
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.1,
"MaxSize": 0.2,
"MeshingDim": 1,
"MinSize": 0.1,
"SelectedEntities": [ 5, 6 ]
},
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.0,
"MaxSize": 2,
"MeshingDim": 2,
"MinSize": 1,
"SelectedEntities": []
}
],
"ThreadNum": 1
}

以下分别为导入的cad模型和最终生成的网格模型。

更多复杂配置请看:二维边界层详细 JSON 示例

Json示例-空间边界层

以三维翼型生成空间边界层混合网格为例。

{
"Distribution": [],
"BoundaryLayer": [
{
"DistributionType": "Separation",
"EntityType": "Face",
"FSPairs":
[
{
"F": 1,
"S": 1
},
{
"F": 8,
"S": 1
},
{
"F": 9,
"S": 1
}
],
"Basic":
{
"GrowthRate":
{
"Constant": 1.2
},
"InitLayerHeight":
{
"Absolute": 0.002
},
"NLayers": 20,
"ProblematicAreasTreatment": "STOP"
},
"SideTreatment":
{
"AutoConnectAngleLimit":true,
"ConnectLimitAngle": 60
}
}
],
"MeshSize": [
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.3,
"MaxSize": 0.2,
"MeshingDim": 1,
"MinSize": 0.1,
"SelectedEntities": [ 1, 2, 3, 4, 17, 18 ]
},
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.1,
"MaxSize": 5,
"MeshingDim": 2,
"MinSize": 1,
"SelectedEntities": [ 2, 3, 4, 5, 6, 7 ]
},
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.1,
"MaxSize": 5,
"MeshingDim": 2,
"MinSize": 0.1,
"SelectedEntities": [ 1, 8, 9 ]
},
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.1,
"MaxSize": 5,
"MeshingDim": 3,
"MinSize": 1,
"SelectedEntities": []
}
],
"ThreadNum": 1
}

以下分别为导入的cad模型和最终生成的网格模型以及剖面图。

Json示例-扫掠网格

注意事项:扫掠体网格只支持单一扫掠体,source源面可以包含多个面,target目标面只能包含单个面。侧面要求是2/3/4边面。

{
"Sweep": [
{
"EntityType": "Face",
"Domain": 1,
"Source": [9],
"Target": [10]
}
],
"MeshSize": [
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.0,
"MaxSize": 1,
"MeshingDim": 3,
"MinSize": 0.1,
"SelectedEntities": []
}
],
"ThreadNum": 1
}

以下分别为导入的cad模型和最终生成的网格模型。

Json示例-复制网格

可以将单个或者多个Entity(边或者面)的网格复制到单个目标Entity上。应用在要求源和目标的网格拓扑完全一致的场景,如周期网格;注意事项:源可以包含多个Entity,目标只能包含一个Entity;源和目标的拓扑结构必须完全一致;当复制面网格时,如果源面和目标面拓扑一致,无须指定边界对应关系"EdgeCorrespondence";如果存在多条边界边对应单一边界边情形,必须正确指定边界边的对应关系。

{
"Copy": [
{
"EntityType": "Face",
"Source": [ 43 ],
"Target": [ 47 ]
}
],
"MeshSize": [
{
"CurvatureFactor": 0,
"GrowthRate": 1.0,
"MaxSize": 3,
"MeshingDim": 2,
"MinSize": 1,
"SelectedEntities": [ 43, 47 ]
}
],
"ThreadNum": 1
}

以Entity为面为例,以下分别为导入的cad模型和复制网格。

Json示例-多种线网格剖分方式

方式一:指定边上单元的数目;
方式二:显式地指定边上每段单元的长度比例;
方式三:按照预设定的函数(Linear/Geometric)分布剖分;

{
"Distribution": [
{
"EntityType": "Edge",
"SelectedEntities": [ 1 ],
"DistributionType": "NumberOfElements",
"NElement": 20
},
{
"EntityType": "Edge",
"SelectedEntities": [ 2 ],
"DistributionType": "Explicit",
"Orientation": "Forward",
"DistributionPara": [ 0.1, 0.3, 0.6, 0.8 ]
},
{
"EntityType": "Edge",
"SelectedEntities": [ 3 ],
"DistributionType": "Linear",
"Orientation": "Forward",
"Symmetric": true,
"NElement": 9,
"EndBeginElementRatio": 5
},
{
"EntityType": "Edge",
"SelectedEntities": [ 4 ],
"DistributionType": "Geometric",
"Orientation": "Forward",
"Symmetric": false,
"NElement": 9,
"EndBeginElementRatio": 5
}
],
"MeshSize": [
{
"CurvatureFactor": 1,
"GrowthRate": 1.0,
"MaxSize": 5,
"MeshingDim": 2,
"MinSize": 0.1,
"SelectedEntities": [1]
}
],
"ThreadNum": 1
}


是否必需 取值范围 默认值 备注
Distribution EntityType Edge
Distribution SelectedEntities 可以为空 Edge tags
Distribution DistributionType NumberOfElements, Explicit, Linear, Geometric 分布类型
Distribution Orientation Forward, Reverse Forward 朝向
Distribution Symmetric true, false false 是否对称分布
Distribution NElement >=1 分段数
Distribution EndBeginElementRatio >0 分布尾段与首段长度之比
Distribution DistributionPara (0.0, 1.0) 严格递增数列

以下为多种线网格剖分结果。

Json示例-高阶网格

支持生成 2 阶高阶网格,导出 VTK 或者 Gmsh 格式文件。

{
"MeshSize": [
{
"CurvatureFactor": 0.2,
"GrowthRate": 1.0,
"MaxSize": 5.0,
"MeshingDim": 3,
"MinSize": 2.0,
"SelectedEntities": []
}
],
"ElementOrder": {
"Order": 2,
"OnGeom":true,
"Serendipity":true
},
"ThreadNum": 1
}


key 是否必需 取值范围 默认值 备注
ElementOrder Order [1,2]
ElementOrder OnGeom [true,false] true 是否贴体
ElementOrder Serendipity [true,false] true true表示高阶点仅存在边上

以下为高阶网格。

导入模型

导入 STEP 文件。

NMAPIModel nmapi;
AMCAX::TopoShape step_shape;
AMCAX::STEP::StepDataTool::Read(step_shape, cadstr);
nmapi.ImportModel({step_shape});
NextMesh 中的模型类
定义 NMAPIModel.hpp:20
AMCAX_API void ImportModel(const std::vector< NMShape > &shapes, const bool replace=true)
加载 AMCAX TopoShape 表示的 CAD 模型。注意:此导入方法不会指定任何实体标签,标签将由内部生成。标签在各维度中从 1 开始计数
static AMCAX_API bool Read(AMCAX::TopoShape &s, std::istream &is)
Read a TopoShape from a stream in STEP format.
形状的基类,包含具有位置和方向信息的基础形状
定义 TopoShape.hpp:15

导入 STL 文件。

const std::string cadstr = "./data/escapement.stl";
NMAPIModel nmapi;
nmapi.ImportSTL(cadstr, Constants::pi / 6.0, Constants::pi/ 3.0, 0.02, 1e-8);
AMCAX_API void ImportSTL(const std::string &filePath, const double dihedralAngle, const double curveAngle, const double mergeSmallPatchRatio=0.0, const double tolerance=0.0)
加载 stl 模型
constexpr double pi
数学常数 Pi,圆的周长与直径之比
定义 Constants.hpp:42

参数说明如下:

选项 取值范围 含义
dihedralAngle [0, Constants::pi] 如果一条边的两个相邻三角面片的二面角大于 dihedralAngle,该边会被检测为特征边
curveAngle [0, Constants::pi] 在相连的特征边(Edge)中,若两个相邻的离散边夹角大于 curveAngle,则特征边(Edge)会从离散边的公共点处拆分为两个子边(Edge)
mergeSmallPatchRatio [0, 1] 假设所有特征边总长为 l,若分割出的 Face 边界周长小于 mergeSmallPatchRatio * l,则会与相邻的面片合并
tolerance [0, inf] 如果两个顶点每个轴的坐标差都小于 tolerance,则两个顶点会合并为一个点

生成/获取网格

导入网格生成控制参数json文件,生成网格。

nlohmann::json paraJ = nlohmann::json::parse(std::ifstream(jsonPath));
nmapi.GenerateMesh(paraJ.dump());
auto meshapi = nmapi.GetMesh();
AMCAX_API void GenerateMesh(const std::string &configJson)
根据配置参数为当前模型生成网格
AMCAX_API NMMesh GetMesh()
获取网格句柄。注意一个模型只有一个唯一的网格

获取网格和几何信息

获取模型实体列表,实体的BRep关系和包围盒,每个实体包含的网格类型与网格单元和顶点,获取网格单元和顶点的全局Id等相关信息。

for (DimType dim : {DimType::D0, DimType::D1, DimType::D2, DimType::D3})
{
std::vector<NMEntity> etVec;
// 获取dim维度下所有实体
nmapi.GetEntities(etVec, dim);
for (auto ent : etVec)
{
auto eTag = nmapi.EntityGetTag(ent);
auto dim = nmapi.EntityGetDim(ent);
// 获取包含指定Entity的dim+1维度实体序列
std::vector<NMEntity> parentVec;
nmapi.GetParentAdjacentEntities(parentVec, ent);
// 获取指定Entity包含的dim-1维度实体序列
std::vector<NMEntity> childVec;
std::vector<Orientation> ories;
nmapi.GetChildAdjacentEntities(childVec, ories, ent);
// 遍历指定Entity包含的Element
auto entElemNum = meshapi.EntityGetElementCount(ent);
for (size_t i = 0; i < entElemNum; i++)
{
auto elem = meshapi.EntityGetElementByIndex(ent, i);
}
// 遍历指定Entity包含的Node
auto entNodeNum = meshapi.EntityGetNodeCount(ent);
for (size_t i = 0; i < entNodeNum; i++)
{
auto node = meshapi.EntityGetNodeByIndex(ent, i);
}
}
}
NMPoint3 pmin, pmax;
// 获取整个模型的boundingbox
nmapi.GetBBox(pmin, pmax);
// 遍历整个Mesh包含的Element
auto elemTotalNum = meshapi.GetElementCount();
for (size_t i = 0; i < elemTotalNum; i++)
{
// 获取第i个Element
auto elem = meshapi.GetElementByIndex(i);
// 获取Element的全局ID,单元类型,包含的Node数
auto elemId = meshapi.ElementGetID(elem);
auto elemType = meshapi.ElementGetType(elem);
auto nodeCount = meshapi.ElementGetNodeCount(elem);
// 遍历指定Element上的所有Node
for (size_t in = 0; in < nodeCount; in++)
{
// 获取Element上的第i个node
auto node = meshapi.ElementGetNode(elem, in);
// 获取Node的全局ID,所属Entity,空间位置,参数
auto nodeId = meshapi.NodeGetID(node);
auto nodeEnt = meshapi.NodeGetEntity(node);
auto nodePos = meshapi.NodeGetPosition(node);
// 仅对边和面内的Node有效
if (nmapi.EntityGetDim(nodeEnt) == DimType::D1 ||
nmapi.EntityGetDim(nodeEnt) == DimType::D2)
double u = meshapi.NodeGetFirstParameter(node);
// 仅对面内的Node有效
if (nmapi.EntityGetDim(nodeEnt) == DimType::D2)
double v = meshapi.NodeGetSecondParameter(node);
}
}
// 遍历整个Mesh包含的Node
auto nodeTotalNum = meshapi.GetNodeCount();
for (size_t i = 0; i < nodeTotalNum; i++)
{
auto node = meshapi.GetNodeByIndex(i);
}
AMCAX_API void GetEntities(std::vector< NMEntity > &ents, const DimType dim)
获取指定维度下的所有实体
AMCAX_API EntTag EntityGetTag(const NMEntity &ent)
获取给定实体的标签
AMCAX_API DimType EntityGetDim(const NMEntity &ent)
获取给定实体的维度
AMCAX_API void GetChildAdjacentEntities(std::vector< NMEntity > &children, std::vector< Orientation > &oris, const NMEntity &ent)
获取包含当前实体的所有 dim-1 维度实体
AMCAX_API void GetParentAdjacentEntities(std::vector< NMEntity > &parents, const NMEntity &ent)
获取包含当前实体的所有 dim+1 维度实体
AMCAX_API void GetBBox(NMPoint3 &pmin, NMPoint3 &pmax, const NMEntity &ent=nullptr)
获取实体的包围盒(当 ent 为空指针时,返回模型的包围盒)
AMCAX::Coord3d NMPoint3
点/向量类型
定义 Macros.hpp:24
DimType
维度类型
定义 Macros.hpp:70

网格加密

支持三角形/四边形/四面体类型网格的贴体/非贴体加密。

meshapi.Refine(false); // false:非贴体加密一次; true:贴体加密一次;

计算网格质量

输出网格单元在指定商软下的质量指标,质量指标包含纵横比,扭曲度,最大最小角度等指标。

double quality = meshapi.ComputeQuality(elem, QualityType::SkewnessByVolume,CommercialSoftware::Fluent);

设置面组

支持面组信息的增删改查。

nmapi.CreatePhysicalSet(DimType::D2, { 1, 2 }, "inlet");
nmapi.CreatePhysicalSet(DimType::D2, { 3, 4 }, "outlet");
nmapi.CreatePhysicalSet(DimType::D2, { 5 }, "symmetry");
std::vector<EntTag> entTags;
nmapi.GetEntitiesInPhysicalSet(entTags, DimType::D2, "outlet");
std::set<std::string> pNames;
nmapi.GetPhysicalSets(pNames, DimType::D2);
AMCAX_API void GetEntitiesInPhysicalSet(std::vector< EntTag > &entTags, const DimType dim, const std::string &pName)
获取物理属性组中的实体
AMCAX_API void GetPhysicalSets(std::set< std::string > &pNames, const DimType dim)
获取指定维度下的所有物理属性组
AMCAX_API void CreatePhysicalSet(const DimType dim, const std::vector< EntTag > &entTags, const std::string &pName)
为指定维度下的实体创建新物理属性组

输出网格文件

输出用户指定格式的网格文件,比如 OBJ/VTK/Fluent MSH/Gmsh MSH。

meshapi.Write("result.vtk", OutFileType::VTK);
meshapi.Write("result.obj", OutFileType::OBJ);
meshapi.Write("result.msh", OutFileType::FLUENT_MSH);
meshapi.WriteGmsh4("result.msh", false, true);

点击这里example01可获得网格剖分示例完整源码,大家根据学习需要自行下载。

附录

下面列出了此示例的完整代码:

#include <fstream>
#include <nlohmann/json.hpp>
#include <set>
using namespace AMCAX::NextMesh;
void GenMesh()
{
// 导入cad文件路径
const std::string cadstr = "./data/qiuwotaojian.stp";
// 网格生成控制参数json文件地址
const std::string jsonPath = "./data/mesh_para-mesh000.json";
// 配置log格式
// 新建NMAPIModel对象
NMAPIModel nmapi;
// 导入cad文件,新建model
AMCAX::TopoShape step_shape;
AMCAX::STEP::StepDataTool::Read(step_shape, cadstr);
nmapi.ImportModel({step_shape});
nlohmann::json paraJ = nlohmann::json::parse(std::ifstream(jsonPath));
// 根据json字符(网格生成控制参数),生成网格
nmapi.GenerateMesh(paraJ.dump());
// 得到生成的网格api
auto meshapi = nmapi.GetMesh();
for (DimType dim : {DimType::D0, DimType::D1, DimType::D2, DimType::D3})
{
std::vector<NMEntity> etVec;
// 获取dim维度下所有实体
nmapi.GetEntities(etVec, dim);
for (auto ent : etVec)
{
auto eTag = nmapi.EntityGetTag(ent);
auto dim = nmapi.EntityGetDim(ent);
// 获取包含指定Entity的dim+1维度实体序列
std::vector<NMEntity> parentVec;
nmapi.GetParentAdjacentEntities(parentVec, ent);
// 获取指定Entity包含的dim-1维度实体序列
std::vector<NMEntity> childVec;
std::vector<Orientation> ories;
nmapi.GetChildAdjacentEntities(childVec, ories, ent);
// 遍历指定Entity包含的Element
auto entElemNum = meshapi.EntityGetElementCount(ent);
for (size_t i = 0; i < entElemNum; i++)
{
auto elem = meshapi.EntityGetElementByIndex(ent, i);
}
// 遍历指定Entity包含的Node
auto entNodeNum = meshapi.EntityGetNodeCount(ent);
for (size_t i = 0; i < entNodeNum; i++)
{
auto node = meshapi.EntityGetNodeByIndex(ent, i);
}
}
}
NMPoint3 pmin, pmax;
// 获取整个模型的boundingbox
nmapi.GetBBox(pmin, pmax);
// 遍历整个Mesh包含的Element
auto elemTotalNum = meshapi.GetElementCount();
for (size_t i = 0; i < elemTotalNum; i++)
{
// 获取第i个Element
auto elem = meshapi.GetElementByIndex(i);
// 获取Element的全局ID,单元类型,包含的Node数
auto elemId = meshapi.ElementGetID(elem);
auto elemType = meshapi.ElementGetType(elem);
auto nodeCount = meshapi.ElementGetNodeCount(elem);
// 遍历指定Element上的所有Node
for (size_t in = 0; in < nodeCount; in++)
{
// 获取Element上的第i个node
auto node = meshapi.ElementGetNode(elem, in);
// 获取Node的全局ID,所属Entity,空间位置,参数
auto nodeId = meshapi.NodeGetID(node);
auto nodeEnt = meshapi.NodeGetEntity(node);
auto nodePos = meshapi.NodeGetPosition(node);
// 仅对边和面内的Node有效
if (nmapi.EntityGetDim(nodeEnt) == DimType::D1 ||
nmapi.EntityGetDim(nodeEnt) == DimType::D2)
double u = meshapi.NodeGetFirstParameter(node);
// 仅对面内的Node有效
if (nmapi.EntityGetDim(nodeEnt) == DimType::D2)
double v = meshapi.NodeGetSecondParameter(node);
}
}
// 遍历整个Mesh包含的Node
auto nodeTotalNum = meshapi.GetNodeCount();
for (size_t i = 0; i < nodeTotalNum; i++)
{
auto node = meshapi.GetNodeByIndex(i);
}
// 设置面组
nmapi.CreatePhysicalSet(DimType::D2, { 1, 2 }, "inlet");
nmapi.CreatePhysicalSet(DimType::D2, { 3, 4 }, "outlet");
nmapi.CreatePhysicalSet(DimType::D2, { 5 }, "symmetry");
std::vector<EntTag> entTags;
nmapi.GetEntitiesInPhysicalSet(entTags, DimType::D2, "outlet");
std::set<std::string> pNames;
nmapi.GetPhysicalSets(pNames, DimType::D2);
// 输出网格文件
meshapi.Write("result.vtk", OutFileType::VTK);
meshapi.Write("result.obj", OutFileType::OBJ);
meshapi.Write("result.msh", OutFileType::FLUENT_MSH);
meshapi.WriteGmsh4("result.msh", false, true);
}
NextMesh 中的模型类
Utility class for operations on STEP shape data structures.