OpenFOAM 编程 | One-Dimensional Transient Heat Conduction
0. 写在前面
本文中将对一维瞬态热传导问题进行数值求解,并基于OpenFOAM类库编写求解器。该问题参考自教科书\(^{[1]}\)示例 8.1。
1. 问题描述
一维瞬态热传导问题控制方程如下
\]
其中,\(\rho c = 1.0\times10^{+7}\ \mathrm{J/m^3\cdot K}\),\(k=10\ \mathrm{W/m\cdot K}\)。
假设等截面直杆长为 \(L=0.02\ \mathrm{m}\),截面为边长 \(0.001\ \mathrm{m}\) 的正方形,全杆初始温度为 \(200 \mathrm{K}\) ,左侧边界条件为\(\nabla T = 0\),右侧边界条件为\(T=0\);杆长方向与 \(x\) 轴平行,此处一维问题不考虑 \(y\) 和 \(z\) 方向。
该问题存在解析解
\]
其中,\(\lambda_n = \frac{\left(2n-1\right)\pi}{2L}\),\(\alpha = \frac{k}{\rho c}\)。
2. 数值解法
对于该物理模型,采用均匀正六面体结构化网格,网格数量为10,相邻网格体心距离为 \(\Delta x = 0.002 \mathrm{m}\),截面面积为\(S = 1\times 10^{-6} \mathrm{m}^2\),网格体积为 \(V_P = 2\times10^{-9} \mathrm{m}^3\),网格示意图如下。

对控制方程进行离散(时间项一阶隐式欧拉格式,固定时间步\(\Delta t = 0.001 \mathrm{s}\)(满足稳定条件)、界面插值采用中心差分格式),可以得到下面线性方程
\]
其对应该问题线性方程组的矩阵形式如下
20.005 & -0.005 & 0 & \cdots & 0 & 0 & 0 \\
-0.005 & 20.005 & -0.005 & \cdots & 0 & 0 & 0\\
\vdots & \vdots & \vdots & \ddots & \vdots & \vdots & \vdots\\
0 & 0 & 0 & \cdots & -0.005 & 20.005 & -0.005 \\
0 & 0 & 0 & \cdots & 0 & -0.005 & 20.015\\
\end{bmatrix}
\begin{bmatrix}
T_0 \\ T_1 \\ \vdots \\ T_8 \\ T_9
\end{bmatrix}
=
\begin{bmatrix}
20 T_0^0 \\ 20 T_1^0 \\ \vdots \\ 20 T_8^0 \\20 T_9^0
\end{bmatrix}
\]
3. OpenFOAM求解
此处,我们把OpenFOAM作为类库使用,可以很快的完成一个求解器,不会涉及过多的底层工作。
3.1 求解器源码
#include "fvCFD.H"
#include <iostream>
int main(int argc, char* argv[])
{
#include "setRootCaseLists.H"
#include "createTime.H" // 构造 runTime 对象
#include "createMesh.H" // 构造 mesh 对象
// 密度 x 热容
dimensionedScalar rhoC("rhoC", dimensionSet(1, -1, -2, 1, 0, 0, 0), scalar(1.e+7));
// 热导率
dimensionedScalar k("k", dimensionSet(1, 1, -3, 1, 0, 0, 0), scalar(10.0));
// 温度场,需要从0文件夹中读取初始值
volScalarField T(IOobject("T", "0", mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);
while ( runTime.loop() )
{
Info << "当前时间 : " << runTime.timeName() << " s" << endl << endl;
// 构造线性方程组
fvScalarMatrix TEqn(fvm::ddt(rhoC, T) == fvm::laplacian(k, T));
// 求解
TEqn.solve();
// 更新边界值
T.correctBoundaryConditions();
if ( runTime.timeIndex() == 1 )
{ // 打印方程组;这段代码放在哪里无所谓,此代码没有在时间步内再次更新 TEqn
Info << "#### UPPER\n" << TEqn.upper() << endl;
Info << "#### DIAG \n" << TEqn.D() << endl;
Info << "#### LOWER\n" << TEqn.lower() << endl;
Info << "#### SOURCE\n" << TEqn.source() << endl; // Right Hand Side
getchar(); // 此处暂停,按回车继续运行...
}
if ( runTime.writeTime() )
{
runTime.write();
}
runTime.printExecutionTime(Info);
}
return 0;
}
3.2 CMakeLists.txt
cmake_minimum_required (VERSION 3.8)
project(OneDimUnsteadyFlow)
# OpenFOAM 安装路径
set( FOAM_PREFIX "/opt/OpenFOAM-v2112" )
# 包含路径
set( FOAM_SRC ${FOAM_PREFIX}/OpenFOAM-v2112/src )
include_directories(
${FOAM_SRC}/atmosphericModels/lnInclude
${FOAM_SRC}/combustionModels/lnInclude
${FOAM_SRC}/conversion/lnInclude
${FOAM_SRC}/dummyThirdParty/lnInclude
${FOAM_SRC}/dynamicFaMesh/lnInclude
${FOAM_SRC}/dynamicFvMesh/lnInclude
${FOAM_SRC}/dynamicMesh/lnInclude
${FOAM_SRC}/engine/lnInclude
${FOAM_SRC}/faOptions/lnInclude
${FOAM_SRC}/fileFormats/lnInclude
${FOAM_SRC}/finiteArea/lnInclude
${FOAM_SRC}/finiteVolume/lnInclude
${FOAM_SRC}/functionObjects/lnInclude
${FOAM_SRC}/fvAgglomerationMethods/lnInclude
${FOAM_SRC}/fvMotionSolver/lnInclude
${FOAM_SRC}/genericPatchFields/lnInclude
${FOAM_SRC}/lagrangian/lnInclude
${FOAM_SRC}/lumpedPointMotion/lnInclude
${FOAM_SRC}/mesh/lnInclude
${FOAM_SRC}/meshTools/lnInclude
${FOAM_SRC}/ODE/lnInclude
${FOAM_SRC}/OpenFOAM/lnInclude
${FOAM_SRC}/optimisation/lnInclude
${FOAM_SRC}/OSspecific/POSIX/lnInclude
${FOAM_SRC}/overset/lnInclude
${FOAM_SRC}/parallel/lnInclude
${FOAM_SRC}/phaseSystemModels/lnInclude
${FOAM_SRC}/Pstream/lnInclude
${FOAM_SRC}/randomProcesses/lnInclude
${FOAM_SRC}/regionFaModels/lnInclude
${FOAM_SRC}/regionModels/lnInclude
${FOAM_SRC}/renumber/lnInclude
${FOAM_SRC}/rigidBodyDynamics/lnInclude
${FOAM_SRC}/rigidBodyMeshMotion/lnInclude
${FOAM_SRC}/sampling/lnInclude
${FOAM_SRC}/semiPermeableBaffle/lnInclude
${FOAM_SRC}/sixDoFRigidBodyMotion/lnInclude
${FOAM_SRC}/sixDoFRigidBodyState/lnInclude
${FOAM_SRC}/surfMesh/lnInclude
${FOAM_SRC}/thermophysicalModels/lnInclude
${FOAM_SRC}/topoChangerFvMesh/lnInclude
${FOAM_SRC}/transportModels/lnInclude
${FOAM_SRC}/TurbulenceModels/lnInclude
${FOAM_SRC}/waveModels/lnInclude
.
)
link_directories(
${FOAM_PREFIX}/ThirdParty-v2112/platforms/linux64Gcc/boost_1_74_0/lib64
${FOAM_PREFIX}/ThirdParty-v2112/platforms/linux64Gcc/fftw-3.3.10/lib64
${FOAM_PREFIX}/ThirdParty-v2112/platforms/linux64Gcc/kahip-3.14/lib64
${FOAM_PREFIX}/ThirdParty-v2112/platforms/linux64GccDPInt32/lib
${FOAM_PREFIX}/ThirdParty-v2112/platforms/linux64GccDPInt32/lib/sys-openmpi
${FOAM_PREFIX}/OpenFOAM-v2112/platforms/linux64GccDPInt32Opt/lib
${FOAM_PREFIX}/OpenFOAM-v2112/platforms/linux64GccDPInt32Opt/lib/dummy
${FOAM_PREFIX}/OpenFOAM-v2112/platforms/linux64GccDPInt32Opt/lib/sys-openmpi
)
set(EXTRA_LIBS dl m)
set(LIBS
Pstream
OpenFOAM
finiteVolume
meshTools
fileFormats
${EXTRA_LIBS}
)
set( CMAKE_CXX_STANDARD 11 )
set( CMAKE_EXE_LINKER_FLAGS "${CMAKE_EXE_LINKER_FLAGS} -Xlinker --no-as-needed -Xlinker --add-needed" )
add_definitions(-Dlinux64 -DWM_ARCH_OPTION=64 -DWM_DP -DWM_LABEL_SIZE=32 -DNoRepository -m64 -fPIC )
# 添加可执行文件
add_executable (${PROJECT_NAME} "main.cpp")
# 链接库
target_link_libraries(${PROJECT_NAME} ${LIBS})
4. 算例设置
4.1 system/blockMeshDict
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object blockMeshDict;
}
scale 1.0; // memter
length 0.02; // 长度
nx 10; // x 方向 网格数量
ny 1;
nz 1;
vertices
(
(0.0 0.0 0.0)
($length 0.0 0.0)
($length 0.001 0.0)
(0.0 0.001 0.0)
(0.0 0.0 0.001)
($length 0.0 0.001)
($length 0.001 0.001)
(0.0 0.001 0.001)
);
edges
(
);
blocks
(
hex (0 1 2 3 4 5 6 7) ($nx $ny $nz) simpleGgrading (1 1 1)
);
boundary
(
left
{
type patch;
faces
(
(0 4 7 3)
);
}
right
{
type patch;
faces
(
(1 2 6 5)
);
}
other
{
type empty;
faces
(
(2 3 7 6)
(4 5 6 7)
(0 1 5 4)
(1 0 3 2)
);
}
);
4.2 system/controDict
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object controlDict;
}
application OneDimUnsteadyFlow;
startFrom startTime;
startTime 0;
stopAt endTime;
endTime 125;
deltaT 0.001;
writeControl adjustableRunTime;
writeInterval 0.1; // 0.1秒为间隔输出数据
purgeWrite 0;
writeFormat ascii;
writePrecision 6;
writeCompression off;
timeFormat general;
timePrecision 6;
runTimeModifiable no;
4.3 system/fvSchemes
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object fvSchemes;
}
ddtSchemes
{
default Euler; // \phi_t = \frac{\phi - \phi^0}{\Delta t}
}
gradSchemes
{
default Gauss linear; // 基于高斯定理的梯度格式
}
divSchemes
{
default Gauss linear; // \phi_f = 0.5(\phi_P + \phi_N)
}
laplacianSchemes
{
default Gauss linear uncorrected; // linear:中心差分格式;uncorrected:不进行非正交性修正
}
本文中,上述ddtSchemes 和 laplacianSchemes 格式为离散方程时所用格式,具体细节已在前边叙述。
4.4 system/fvSolution
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object fvSolution;
}
solvers
{
T
{
solver GAMG;
smoother GaussSeidel;
tolerance 1e-08;
relTol 0.0;
}
}
4.5 0/T
FoamFile
{
version 2.0;
format ascii;
arch "LSB;label=32;scalar=64";
class volScalarField;
location "0";
object T;
}
dimensions [0 0 0 1 0 0 0];
internalField uniform 200;
boundaryField
{
left
{
type zeroGradient; // 零梯度
}
right
{
type fixedValue; // 固定值
value uniform 0;
}
other
{
type empty;
}
}
5. 求解计算
文件结构如下所示。
.
├── build // build 目录,用于编译代码
├── CMakeLists.txt // 项目管理,内容见3.2节
├── main.cpp // 求解器源代码,内容见3.1节
└── OneDimUnsteadyFlowCase // 算例所在目录 ***
├── 0 // 0 文件夹,保存初始条件
│ └── T // 本示例中只有温度场 T,故此处只有 T 文件,内容见 4.5 节
├── OneDimUnsteadyFlowCase.foam // 算例目录名称+foam扩展名,空文件,仅作ParaView加载结果使用
└── system // system 目录
├── blockMeshDict // blockMesh字典文件,内容见 4.1 节
├── controlDict // 求解器运行控制字典文件,内容见 4.2 节
├── fvSchemes // 有限体积数值格式字典文件,内容见 4.3 节
└── fvSolution // 求解器参数设置字典文件,内容见 4.4 节
5.1 编译求解器

主要命令解释:
$ cd build/ # 从当前目录切换到路径 ./build
$ cmake .. # 执行 CMake 生成构建文件(当前生成的是MakefIle)
$ make # 执行 Make,编译代码
5.2 运行求解器

主要命令解释:
$ cd OneDimUnsteadyFlowCase/ # 从当前目录切换到算例目录
$ blockMesh > log.blockMesh # 运行 blockMesh 画网格,并将标准输出重定向到 log.blockMesh
$ ../build/OneDimUnsteadyFlow # 运行求解器,注意求解器的相对路径
另外,我们也可以看到求解器打印的线性方程组与数值解法中所描述的是一致的。
6. 后处理
我们对比第 \(40 \ \mathrm{s}\) 时数值结果与解析结果
云图:

曲线图:

参考文献
[1] H. Versteeg , W. Malalasekera. Introduction to Computational Fluid Dynamics, An: The Finite Volume Method 2nd Edition[M]. Pearson. 2007
OpenFOAM 编程 | One-Dimensional Transient Heat Conduction的更多相关文章
- OpenFOAM 编程 | 求解捕食者与被捕食者模型(predator-prey model)问题(ODEs)
0. 写在前面 本文问题参考自文献 \(^{[1]}\) 第一章例 6,并假设了一些条件,基于 OpenFOAM-v2206 编写程序数值上求解该问题.笔者之前也写过基于 OpenFOAM 求解偏分方 ...
- CUDA Samples: heat conduction(模拟热传导)
以下CUDA sample是分别用C++和CUDA实现的模拟热传导生成的图像,并对其中使用到的CUDA函数进行了解说,code参考了<GPU高性能编程CUDA实战>一书的第七章,各个文件内 ...
- OpenFOAM编程 | Hello OpenFOAM
写在前面 OpenFOAM 是一个非常好用的开源程序包,笔者一直在研究和使用,其编程语言是笔者非常喜欢使用的 C++.但是笔者不是很喜欢 OpenFOAM 自己的构建工具 wmake,更倾向于使用 C ...
- LibTorch | 使用神经网络求解一维稳态对流扩散方程
0. 写在前面 本文将使用基于LibTorch(PyTorch C++接口)的神经网络求解器,对一维稳态对流扩散方程进行求解.研究问题参考自教科书\(^{[1]}\)示例 8.3. 目录 0. 写在前 ...
- OpenFOAM 学习路线 【转载】
"Two weeks of playing with a CFD code will save you one afternoon of reading" 什么是OpenFOAM( ...
- Java基础面试题(Hibernate)
Hibernate是一个什么样的框架? Hibernate是一个开放源代码的对象关系映射框架,它对JDBC进行了非常轻量级的对象封装,它将POJO与数据库表建立映射关系,是一个全自动的orm框架,hi ...
- Hibernate入门(3)- 持久对象的生命周期介绍
在hibernate中对象有三种状态:瞬时态(Transient). 持久态(Persistent).脱管态或游离态(Detached).处于持久态的对象也称为PO(Persistence Objec ...
- 【Java编程】volatile和transient关键字的理解
理解volatile volatile是java提供的一种轻量级的同步机制,所谓轻量级,是相对于synchronized(重量级锁,开销比较大)而言的. 根据java虚拟机的内存模型,我们知道 ...
- OpenFOAM Tutorial Standard Solvers【转载】
转载自:http://www.cnblogs.com/fortran/articles/1996927.html boundaryFoam Steady-state solver for 1D tur ...
随机推荐
- UiPath选择器之页面选择器的介绍和使用
一.页面选择器的介绍 某些软件程序的布局和属性节点具有易变的值,例如某些Web应用程序.UiPath Studio无法预测这些变化,因此,您可能必须手动生成一些选择器. 每个属性都有一个分配的值.选择 ...
- sql-DDL-约束
约束 对表中的数据进行限定,保证数据的正确性.有效性和完整性. 6个约束 1. 主键约束Primary Key: 唯一,不能为null -- 主键约束.和唯一约束不能同时设置 1. 含义:非空且唯一 ...
- 20行python代码,轻松获取各路小说,非常简单
哔哔两句 作为现代青年,我相信应该没几个没看过小说的吧,嘿嘿~ 一般来说咱们书荒的时候怎么办?自然是去起某点排行榜先找到小说名字,然后再找度娘一搜,哎 ,笔趣阁就出来答案了,美滋滋~但是那多麻烦,咱们 ...
- NC16430 [NOIP2016]蚯蚓
NC16430 [NOIP2016]蚯蚓 题目 题目描述 本题中,我们将用符号 \(\lfloor c \rfloor\) 表示对 c 向下取整,例如:\(\lfloor 3.0 \rfloor = ...
- 详解HashMap源码解析(下)
上文详解HashMap源码解析(上)介绍了HashMap整体介绍了一下数据结构,主要属性字段,获取数组的索引下标,以及几个构造方法.本文重点讲解元素的添加.查找.扩容等主要方法. 添加元素 put(K ...
- 重磅硬核 | 一文聊透对象在 JVM 中的内存布局,以及内存对齐和压缩指针的原理及应用
欢迎关注公众号:bin的技术小屋 大家好,我是bin,又到了每周我们见面的时刻了,我的公众号在1月10号那天发布了第一篇文章<从内核角度看IO模型的演变>,在这篇文章中我们通过图解的方式以 ...
- freeswitch拨打分机号源代码跟踪
概述 freeswitch是一款非常好用的开源VOIP软交换平台. 之前我们有介绍过使用fs拨打分机号的方法,其中代码流程是比较复杂的,所以单独开一章介绍. fs拨打分机号,是使用send_dtmf接 ...
- Unbuntu VS Code启动时闪退暂时的解决方法
背景: 刚刚试着更新了操作系统,没更新成功,在下载系统更新的时候brave浏览器消失了,wps消失了,搜狗拼音输入法消失了.更新时,卡在Kernal Offset上,然后长按电源键再重启就好了.但是v ...
- 一张图进阶 RocketMQ - 通信机制
前 言 三此君看了好几本书,看了很多遍源码整理的 一张图进阶 RocketMQ 图片,关于 RocketMQ 你只需要记住这张图!觉得不错的话,记得点赞关注哦. [重要]视频在 B 站同步更新,欢迎围 ...
- S32K148-CAN收发
最近在搞一个转换板开发,大概意思把CAN信号转换成SPI信号,方案有两种:1)通过硬件电路直接把信号的bit位一位一位移给两个集成芯片:2)通过MCU接收CAN信号,再把信号变量转换成SPI信号发送给 ...