物体から流体への熱伝達(自然対流冷却)

  • シミュレーション名 :: heatTransfer__box_in_water_XY2D

シミュレーション体系

  • 基本方程式は、 Navier-Stokes方程式標準k-ε方程式熱輸送方程式

    • 流体は、 非圧縮 , 標準k-εモデル を取扱い、対象は空気とする.

  • 計算対象は、2次元の閉じた立方体の箱の中にバルク鉄塊をおいた際の熱伝達をみる.

    • 2次元 長さ 200 [mm], 幅 200 [mm]、高さ 200 [mm]、厚み 10 [mm] の中空の箱の中央下部に、長さ 50 [mm] 、幅 50 [mm]、高さ 50 [mm] の鉄塊(SS400)を置く.

    • 標準k-εモデル を用いて、壁境界は 垂直方向流速にゼロを課した (Normal-Tangential).

    • 全領域の初期条件は、流速は0、k=0.00457, epsilon=1e-4 とおいている.

物性条件

Materials Settings

Target

Parameters

Value

Unit

Description

空気

Density

1.2e+0

kg/m3

密度

viscosity

1.0e-5

Pa.s

粘度

KE SigmaK

1.0

KE SigmaE

1.3

KE C1

1.44

係数C1

KE C2

1.92

係数C1

KE Cmu

0.09

係数Cmu

KE Clip

1.0e-6

Viscosity Model

K-Epsilon

使用するモデル

Heat Conductivity

0.0257

W/m.K

熱伝導度

Heat Capacity

1.005e+3

J/kg.K

比熱

Density

1.166e+0

kg/m3

密度

SS400

Heat Conductivity

51.6

W/m.K

熱伝導率

Heat Capacity

473.0

J/kg.K

比熱

Density

7.85e+3

kg/m3

密度

メッシュ

  • メッシュ生成スクリプト ( mesh.py )

mesh.py ( heatTransfer__box_in_water_XY2D )
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()
    gmsh.model.occ.removeAllDuplicates()
    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

geometry.conf ( heatTransfer__box_in_water_XY2D )
# <define> @cntnr.Lx      = 200e-3
# <define> @cntnr.Ly      = 200e-3
# <define> @cntnr.thick   =  10e-3
# <define> @cntnr.x0      = (-0.5) * @cntnr.Lx
# <define> @cntnr.y0      = 0.0
# <define> @cntnr.dx	  = @cntnr.Lx
# <define> @cntnr.dy	  = @cntnr.Ly

# <define> @fluid.Lx      = @cntnr.Lx - 2.0 * @cntnr.thick
# <define> @fluid.Ly      = @cntnr.Ly - 2.0 * @cntnr.thick
# <define> @fluid.x0      = (-0.5) * @fluid.Lx
# <define> @fluid.y0      = @cntnr.thick
# <define> @fluid.dx	  = @fluid.Lx
# <define> @fluid.dy	  = @fluid.Ly

# <define> @metal.Lx      = 50e-3
# <define> @metal.Ly      = 50e-3
# <define> @metal.x0      = (-0.5) * @metal.Lx
# <define> @metal.y0      = @cntnr.thick
# <define> @metal.dx	  = @metal.Lx
# <define> @metal.dy	  = @metal.Ly


# <names> key   geometry_type   x0          y0          dx	     dy
metal      	quad		@metal.x0   @metal.y0	@metal.dx    @metal.dy
fluid.00      	quad		@fluid.x0   @fluid.y0	@fluid.dx    @fluid.dy
cntnr.00  	quad		@cntnr.x0   @cntnr.y0	@cntnr.dx    @cntnr.dy

# <names> key   boolean_type    targetKeys     toolKeys	  removeObject removeTool
container 	cut		[container.00] [fluid.00] True	       False
fluid		cut		[fluid.00]     [metal] 	  True	       False
  • boundary.json

boundary.json ( heatTransfer__box_in_water_XY2D )
{
    "metal" : [ [2,1] ],
    "fluid": [ [2,2] ],
    "container": [ [2,3] ],
    "metal_fluid": [ [1,5], [1,6], [1,7] ], 
    "metal_container": [ [1,13] ], 
    "fluid_container": [ [1,1], [1,2], [1,3], [1,4], [1,8] ], 
    "container_outside": [ [1,9], [1,10], [1,11], [1,12] ]
}
  • phys.conf

phys.conf ( heatTransfer__box_in_water_XY2D )
# <names> key	   type		dimtags_keys			physNum
metal	  	   surf		[metal]				201
fluid		   surf		[fluid]				202
container	   surf		[container]			203
metal_fluid	   line		[metal_fluid]			101
metal_container	   line		[metal_container]		102
fluid_container	   line		[fluid_container]		103
container_outside  line		[container_outside]		104
  • mesh.conf

mesh.conf ( heatTransfer__box_in_water_XY2D )
# <names> key	    physNum	  meshType    resolution1	resolution2	evaluation
metal	  	    201	  	  constant    0.01		-		-
fluid	  	    202	  	  constant    0.01		-		-
container	    203	  	  constant    0.01		-		-
metal_fluid  	    101	  	  constant    0.0		-		-
metal_container	    102	  	  constant    0.0		-		-
fluid_container	    103	  	  constant    0.0		-		-
container_outside   104	  	  constant    0.0		-		-
  • 生成したメッシュを次に示す.

../../../_images/mesh3.png

elmer シミュレーション設定

  • elmer シミュレーション設定ファイルを以下に示す.

steady_ke.sif ( heatTransfer__box_in_water_XY2D )
!! ========================================================= !!
!! ===  time-evolution of natural cooling by Boussinesq  === !!
!! ========================================================= !!

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                        = "Cartesian"
  Coordinate Mapping(3)                    = 1 2 3

  Simulation Type                          = "Transient"
  TimeStepping Method                      = BDF
  BDF Order                                = 2
  Timestep sizes(1)                        = 1e-1
  Timestep Intervals(1)                    = 1000

  Steady State Max Iterations              = 200
  Post File                                = heatTransfer.vtu
End


Constants
  Stefan Boltzmann                         = 5.6703e-8    !!  -- [ W / m^2 K^4 ] --  !!
  gravity(4)                               = 0 -1 0 9.82  !!  -- [ m/s^2 ] -- !!
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


Solver 3
  Equation                                 = "HeatEquations"
  Procedure                                = "HeatSolve" "HeatSolver"
  Variable                                 = "Temperature"
  Exec Solver                              = "Always"

  Stabilize                                = True
  Optimize Bandwidth                       = True

  Linear System Solver                     = Iterative
  Linear System Iterative Method           = BiCGStab
  Linear System Max Iterations             = 500
  Linear System Convergence Tolerance      = 1.0e-8
  Linear System Preconditioning            = ILU0
  Linear System Precondition Recompute     = 1

  Steady State Convergence Tolerance       = 1.0e-5
  Nonlinear System Convergence Tolerance   = 1.0e-4
  Nonlinear System Max Iterations          = 1
  Nonlinear System Newton After Iterations = 10000
  Nonlinear System Newton After Tolerance  = 0.0
  Nonlinear System Relaxation Factor       = 0.5

  Steady State Convergence Tolerance       = 1.0e-5
End



Equation 1
  Name                                     = "metals"
  Active Solvers(1)                        = 3
  Convection                               = Computed
End

Equation 2
  Name                                     = "fluid"
  Active Solvers(3)                        = 2 1 3
  Convection                               = Computed
End


!! ------------------------------------------------- !!
!! --- [3] body & materials                      --- !!
!! ------------------------------------------------- !!

Body 1
  Name                                     = "metal"
  Target Bodies(1)                         = $metal
  Equation                                 = 1
  Material                                 = 1
  Initial Condition                        = 1
End

Body 2
  Name                                     = "fluid"
  Target Bodies(1)                         = $fluid
  Equation                                 = 2
  Material                                 = 2
  Initial Condition                        = 2
  Body Force                               = 1
End

Body 3
  Name                                     = "container"
  Target Bodies(1)                         = $container
  Equation                                 = 1
  Material                                 = 1
  Initial Condition                        = 3
End

Material 1
  Name                                     = "SS400"
  Heat Conductivity                        = 51.6       !! W/m.K
  Heat Capacity                            = 473.0      !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Density                                  = 7.85e+3    !! kg/m3
End


Material 2
  Name                                     = "Air"
  !! fluid property !!
  Viscosity                                = 1.0e-5     !! Pa.s
  Density                                  = 1.166e+0   !! kg/m3
  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
  !! heat property !!
  Heat Conductivity                        = 0.0257     !! W/m.K
  Heat Capacity                            = 1.005e+3   !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Heat Expansion Coefficient               = 3.665e-3   !! 1e-6
End

!! ------------------------------------------------- !!
!! --- [4] initial boundary condition            --- !!
!! ------------------------------------------------- !!

$ T_metal     = 473.15
$ T_fluid     = 293.15
$ T_container = 293.15

Initial Condition 1
  Name                                     = "metal.initial"
  Temperature                              = $T_metal
End

Initial Condition 2
  Name                                     = "fluid.initial"
  Velocity 1                               = 0
  Velocity 2                               = 0
  Temperature                              = $T_fluid

  Kinetic Energy                           = 0.00457
  Kinetic Dissipation                      = 1.0e-4
End

Initial Condition 3
  Name                                     = "container.initial"
  Temperature                              = $T_container
End


Boundary Condition 1
  Name                                     = "fluid.boundary"
  Target Boundaries(2)                     = 102 103
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
  !! Noslip Wall BC                           = True
End


!! ------------------------------------------------- !!
!! --- [6] others                                --- !!
!! ------------------------------------------------- !!
Body Force 1
  Name                                     = "BoussinesqApprox"
  Boussinesq                               = Logical True
End

自然対流冷却のシミュレーション結果

温度の時間発展

../../../_images/temperature_anime.gif

結果について

  • 鉄塊内部の温度差が顕著.

  • 容器側へ、熱伝導で大部分の熱が逃げている可能性が高い.

  • 容器なしでシミュレーションして比較検討してみる.

容器なしの場合との比較

メッシュ

  • メッシュ生成スクリプト ( mesh.py )

mesh.py ( heatTransfer__box_in_water_XY2D )
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/woContainer/geometry.conf"
    bdrFile = "dat/woContainer/boundary.json"
    dimtags = dg2.define__geometry_2d( inpFile=inpFile, dimtags=dimtags )
    dimtags = ldt.load__dimtags( dimtags=dimtags, inpFile=bdrFile )
    gmsh.model.occ.synchronize()
    gmsh.model.occ.removeAllDuplicates()
    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/woContainer/mesh.conf"
        physFile = "dat/woContainer/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

geometry.conf ( heatTransfer__box_in_water_XY2D )
# <define> @cntnr.Lx      = 200e-3
# <define> @cntnr.Ly      = 200e-3
# <define> @cntnr.thick   =  10e-3
# <define> @cntnr.x0      = (-0.5) * @cntnr.Lx
# <define> @cntnr.y0      = 0.0
# <define> @cntnr.dx	  = @cntnr.Lx
# <define> @cntnr.dy	  = @cntnr.Ly

# <define> @fluid.Lx      = @cntnr.Lx - 2.0 * @cntnr.thick
# <define> @fluid.Ly      = @cntnr.Ly - 2.0 * @cntnr.thick
# <define> @fluid.x0      = (-0.5) * @fluid.Lx
# <define> @fluid.y0      = @cntnr.thick
# <define> @fluid.dx	  = @fluid.Lx
# <define> @fluid.dy	  = @fluid.Ly

# <define> @metal.Lx      = 50e-3
# <define> @metal.Ly      = 50e-3
# <define> @metal.x0      = (-0.5) * @metal.Lx
# <define> @metal.y0      = @cntnr.thick
# <define> @metal.dx	  = @metal.Lx
# <define> @metal.dy	  = @metal.Ly


# <names> key   geometry_type   x0          y0          dx	     dy
metal      	quad		@metal.x0   @metal.y0	@metal.dx    @metal.dy
fluid.00      	quad		@fluid.x0   @fluid.y0	@fluid.dx    @fluid.dy

# <names> key   boolean_type    targetKeys     toolKeys	  removeObject removeTool
fluid		cut		[fluid.00]     [metal] 	  True	       False
  • boundary.json

boundary.json ( heatTransfer__box_in_water_XY2D )
{
    "metal" : [ [2,1] ],
    "fluid": [ [2,2] ],
    "metal_fluid": [ [1,5], [1,6], [1,7] ], 
    "fluid_outside": [ [1,1], [1,2], [1,3], [1,4], [1,8] ], 
    "metal_outside": [ [1,9] ]
}
  • phys.conf

phys.conf ( heatTransfer__box_in_water_XY2D )
# <names> key	   type		dimtags_keys			physNum
metal	  	   surf		[metal]				201
fluid		   surf		[fluid]				202
metal_fluid	   line		[metal_fluid]			101
fluid_outside  	   line		[fluid_outside]			102
metal_outside  	   line		[metal_outside]			103
  • mesh.conf

mesh.conf ( heatTransfer__box_in_water_XY2D )
# <names> key	    physNum	  meshType    resolution1	resolution2	evaluation
metal	  	    201	  	  constant    0.01		-		-
fluid	  	    202	  	  constant    0.01		-		-
metal_fluid  	    101	  	  constant    0.0		-		-
fluid_outside       102	  	  constant    0.0		-		-
metal_outside	    103		  constant    0.0		-		-

容器なしの場合のシミュレーション設定

  • elmer シミュレーション設定ファイルを以下に示す.

steady_ke.sif ( heatTransfer__box_in_water_XY2D )
!! ========================================================= !!
!! ===  time-evolution of natural cooling by Boussinesq  === !!
!! ========================================================= !!

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                        = "Cartesian"
  Coordinate Mapping(3)                    = 1 2 3

  Simulation Type                          = "Transient"
  TimeStepping Method                      = BDF
  BDF Order                                = 2
  Timestep sizes(1)                        = 5e-2
  Timestep Intervals(1)                    = 500

  Steady State Max Iterations              = 200
  Post File                                = steady_ke.vtu
End

Constants
  Stefan Boltzmann                         = 5.6703e-8    !!  -- [ W / m^2 K^4 ] --  !!
  gravity(4)                               = 0 -1 0 9.82  !!  -- [ m/s^2 ] -- !!
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


Solver 3
  Equation                                 = "HeatEquations"
  Procedure                                = "HeatSolve" "HeatSolver"
  Variable                                 = "Temperature"
  Exec Solver                              = "Always"

  Stabilize                                = True
  Optimize Bandwidth                       = True

  Linear System Solver                     = Iterative
  Linear System Iterative Method           = BiCGStab
  Linear System Max Iterations             = 500
  Linear System Convergence Tolerance      = 1.0e-8
  Linear System Preconditioning            = ILU0
  Linear System Precondition Recompute     = 1

  Steady State Convergence Tolerance       = 1.0e-5
  Nonlinear System Convergence Tolerance   = 1.0e-4
  Nonlinear System Max Iterations          = 1
  Nonlinear System Newton After Iterations = 10000
  Nonlinear System Newton After Tolerance  = 0.0
  Nonlinear System Relaxation Factor       = 0.5

  Steady State Convergence Tolerance       = 1.0e-5
End



Equation 1
  Name                                     = "metals"
  Active Solvers(1)                        = 3
  Convection                               = Computed
End

Equation 2
  Name                                     = "fluid"
  Active Solvers(3)                        = 2 1 3
  Convection                               = Computed
End


!! ------------------------------------------------- !!
!! --- [3] body & materials                      --- !!
!! ------------------------------------------------- !!

Body 1
  Name                                     = "metal"
  Target Bodies(1)                         = $metal
  Equation                                 = 1
  Material                                 = 1
  Initial Condition                        = 1
End

Body 2
  Name                                     = "fluid"
  Target Bodies(1)                         = $fluid
  Equation                                 = 2
  Material                                 = 2
  Initial Condition                        = 2
  Body Force                               = 1
End

Material 1
  Name                                     = "SS400"
  Heat Conductivity                        = 51.6       !! W/m.K
  Heat Capacity                            = 473.0      !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Density                                  = 7.85e+3    !! kg/m3
End


Material 2
  Name                                     = "Air"
  !! fluid property !!
  Viscosity                                = 1.0e-5     !! Pa.s
  Density                                  = 1.166e+0   !! kg/m3
  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
  !! heat property !!
  Heat Conductivity                        = 0.0257     !! W/m.K
  Heat Capacity                            = 1.005e+3   !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Heat Expansion Coefficient               = 3.665e-3   !! 1e-6
End

!! ------------------------------------------------- !!
!! --- [4] initial boundary condition            --- !!
!! ------------------------------------------------- !!

$ T_metal     = 473.15
$ T_fluid     = 293.15

Initial Condition 1
  Name                                     = "metal.initial"
  Temperature                              = $T_metal
End

Initial Condition 2
  Name                                     = "fluid.initial"
  Velocity 1                               = 0
  Velocity 2                               = 0
  Temperature                              = $T_fluid

  Kinetic Energy                           = 0.00457
  Kinetic Dissipation                      = 1.0e-4
End


Boundary Condition 1
  Name                                     = "fluid.boundary"
  Target Boundaries(2)                     = 101 102
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
  !! Noslip Wall BC                           = True
End


!! ------------------------------------------------- !!
!! --- [6] others                                --- !!
!! ------------------------------------------------- !!
Body Force 1
  Name                                     = "BoussinesqApprox"
  Boussinesq                               = Logical True
End

容器の有無による温度分布の差異 ( t=25 (s) 経過後 )

../../../_images/results_T.png

比較した結果について

  • 明らかに容器を伝って熱伝導により、熱が逃げているのがわかる.