標準 k-ε モデルを用いたバックステップ流れ¶
シミュレーション名 :: flow__backStep_keSolver_XY2D
シミュレーション体系¶
基本方程式は、 Navier-Stokes方程式 と 標準k-ε方程式
流体は、 非圧縮 , 標準k-εモデル を取扱い、対象は空気とする.
計算対象は、いわゆるバックステップ流れと呼ばれるものである.
2次元 長さ 13 [m], 幅 1 [m] の流路から、長さ 13 [m] 、幅 2 [m] の流路への切り替わる.
流入境界は、ポワズイユ流れとし、流出境界は uy=0 とする.
標準k-εモデル を用いて、壁境界は 滑りなし壁条件を課した.
全領域の初期条件は、流速は0、k=0.00457, epsilon=1e-4 とおいている.
メッシュ¶
メッシュ生成スクリプト ( mesh.py )
import numpy as np
import os, sys
import gmsh
import nkGmshRoutines.define__geometry_2d as dg2
import nkGmshRoutines.load__dimtags as ldt
# ========================================================= #
# === 実行部 === #
# ========================================================= #
if ( __name__=="__main__" ):
# ------------------------------------------------- #
# --- [1] initialization of the gmsh --- #
# ------------------------------------------------- #
gmsh.initialize()
gmsh.option.setNumber( "General.Terminal", 1 )
gmsh.option.setNumber( "Mesh.Algorithm" , 5 )
gmsh.option.setNumber( "Mesh.Algorithm3D", 4 )
gmsh.option.setNumber( "Mesh.SubdivisionAlgorithm", 0 )
gmsh.model.add( "model" )
# ------------------------------------------------- #
# --- [2] Modeling --- #
# ------------------------------------------------- #
dimtags = {}
inpFile = "dat/geometry.conf"
bdrFile = "dat/boundary.json"
dimtags = dg2.define__geometry_2d( inpFile=inpFile, dimtags=dimtags )
dimtags = ldt.load__dimtags( dimtags=dimtags, inpFile=bdrFile )
gmsh.model.occ.synchronize()
# ------------------------------------------------- #
# --- [3] Mesh settings --- #
# ------------------------------------------------- #
mesh_from_config = True # from nkGMshRoutines/test/mesh.conf, phys.conf
uniform_size = 0.05
if ( mesh_from_config ):
meshFile = "dat/mesh.conf"
physFile = "dat/phys.conf"
import nkGmshRoutines.assign__meshsize as ams
meshes = ams.assign__meshsize( meshFile=meshFile, physFile=physFile, \
dimtags=dimtags, target="surf" )
else:
import nkGmshRoutines.assign__meshsize as ams
meshes = ams.assign__meshsize( uniform=uniform_size, dimtags=dimtags )
# ------------------------------------------------- #
# --- [4] post process --- #
# ------------------------------------------------- #
gmsh.model.occ.synchronize()
gmsh.model.mesh.generate(2)
gmsh.write( "msh/model.msh" )
gmsh.finalize()
geometry.conf
# <define> @pt1 = [ 0,2,0]
# <define> @pt2 = [26,2,0]
# <define> @pt3 = [26,0,0]
# <define> @pt4 = [13,0,0]
# <define> @pt5 = [13,1,0]
# <define> @pt6 = [ 0,1,0]
# <names> key geometry_type point1 point2 point3 point4 point5 point6
fluid point_surf @pt1 @pt2 @pt3 @pt4 @pt5 @pt6
boundary.json
{
"inlet" : [ [1,6] ],
"outlet": [ [1,2] ],
"wall_x": [ [1,1], [1,3], [1,5] ],
"wall_y": [ [1,4] ]
}
phys.conf
# <names> key type dimtags_keys physNum
fluid surf [fluid] 201
inlet line [inlet] 101
outlet line [outlet] 102
wall_x line [wall_x] 103
wall_y line [wall_y] 104
mesh.conf
# <names> key physNum meshType resolution1 resolution2 evaluation
fluid 201 constant 0.05 - -
inlet 101 constant 0.0 - -
outlet 102 constant 0.0 - -
wall_x 103 constant 0.0 - -
wall_y 104 constant 0.0 - -
生成したメッシュを次に示す.
elmer シミュレーション設定¶
elmer シミュレーション設定ファイルを以下に示す.
!! ========================================================= !!
!! === steady state of back-step flow === !!
!! ========================================================= !!
include "./msh/model/mesh.names"
Header
CHECK KEYWORDS Warn
Mesh DB "." "msh/model"
Include Path ""
Results Directory "out/"
End
!! ------------------------------------------------- !!
!! --- [1] Simulation --- !!
!! ------------------------------------------------- !!
Simulation
Max Output Level = 3
Coordinate System = string "Cartesian"
Coordinate Mapping(3) = 1 2 3
Simulation Type = "Steady State"
Steady State Max Iterations = 200
Post File = steady_ke.vtu
End
!! ------------------------------------------------- !!
!! --- [2] Solver & Equation --- !!
!! ------------------------------------------------- !!
Solver 1
Equation = k-epsilon
Procedure = "KESolver" "KESolver"
Stabilize = True
Linear System Solver = Direct
Linear System Direct Method = String "umfpack"
!! Linear System Iterative Method = BiCGStab
Linear System Max Iterations = 10000
Linear System Preconditioning = ILUT
Linear System Convergence Tolerance = 1.0e-5
Nonlinear System Max Iterations = 1
Nonlinear System Convergence Tolerance = 1.0e-5
Nonlinear System Relaxation Factor = 0.5
Nonlinear System Newton After Tolerance = 0.0
Nonlinear System Newton After Iterations = 10000
Steady State Convergence Tolerance = 1.0e-5
End
Solver 2
Equation = Navier-Stokes
Procedure = "FlowSolve" "FlowSolver"
Stabilize = True
Linear System Solver = Direct
Linear System Direct Method = String "umfpack"
!! Linear System Iterative Method = BiCGStab
Linear System Max Iterations = 10000
Linear System Convergence Tolerance = 1.0e-5
Linear System Preconditioning = ILUT
Nonlinear System Max Iterations = 1
Nonlinear System Convergence Tolerance = 1.0e-5
Nonlinear System Relaxation Factor = 0.5
Nonlinear System Newton After Tolerance = 0.0
Nonlinear System Newton After Iterations = 10000
Steady State Convergence Tolerance = 1.0e-5
End
Equation 1
Name = "Fluid"
Active Solvers(2) = 2 1
End
!! ------------------------------------------------- !!
!! --- [3] body & materials --- !!
!! ------------------------------------------------- !!
Body 1
Name = "Fluid"
Target Bodies(1) = $ fluid
Equation = 1
Material = 1
Initial Condition = 1
End
Material 1
Name = "Air"
Viscosity = 1.0e-5
Density = 1.2e0
KE SigmaK = 1.00
KE SigmaE = 1.30
KE C1 = 1.44
KE C2 = 1.92
KE Cmu = 0.09
KE Clip = Real 1.0e-6
Viscosity Model = K-Epsilon
End
!! ------------------------------------------------- !!
!! --- [4] initial boundary condition --- !!
!! ------------------------------------------------- !!
Initial Condition 1
Velocity 1 = 0
Velocity 2 = 0
Kinetic Energy = 0.00457
Kinetic Dissipation = 1.0e-4
End
Boundary Condition 1
Name = "inlet"
Target Boundaries(1) = $ inlet
Velocity 1 = Variable Coordinate 2
Real MATC "6*(tx-1)*(2-tx)"
Velocity 2 = 0
Kinetic Energy = Real 0.00457
Kinetic Dissipation = Real 1.0e-4
End
Boundary Condition 2
Name = "outlet"
Target Boundaries(1) = $ outlet
Velocity 2 = 0.0
End
Boundary Condition 3
Name = "wall_x"
Target Boundaries(1) = $ wall_x
Noslip Wall BC = True
End
Boundary Condition 4
Name = "wall_y"
Target Boundaries(1) = $ wall_y
Noslip Wall BC = True
End
!! Solver 2 :: Reference Norm = Real 0.23583018
!! Solver 2 :: Reference Norm Tolerance = Real 1.0e-3
!! RUN
バックステップ流れのシミュレーション¶
バックステップ流れの定常解 ( 標準 k-ε モデル )¶
Steady State 定常解を求めると、80回程度の反復で閾値( 1e-5 )以下に収束した.
反復法を用いた場合、収束性が問題となり、最適な問題設定(閾値、前処理、その他、)ができず、行き詰まった.
N-Sソルバ単体だとやはり、計算が爆発する... ( k-εモデルの必要性? : uuの項が高波数モードをだす? )
線形解法に直接法 ( Direct, umfpack ) を用いると、収束した.
N-S全般で、収束性が問題となっているので、N-Sには直接法を用いるのがよいかもしれない.
結果について¶
配管中の速度はほぼ流入設定と同様のポワズイユ流れとなっている.
段差部分において、流れの剥離が生じる.
圧力は、細管部分で主に損失し、太管部分は変化に乏しい.
kinetic energy ( 渦のエネルギー )は、配管壁付近(境界層付近)に集中する.
dissipation も同様.散逸は、壁付近で生じる.