水中に沈めた鉄塊からの自然対流冷却

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

シミュレーション体系

  • 基本方程式は、 熱輸送方程式( Heat Equation )Navier-Stokes方程式

    • 流体は、 非圧縮 , Boussinesq近似

    • 熱輸送は、 熱伝導、熱対流、輻射を取り扱う.

  • シミュレーション対象は直方体 ( 50 [cm] x 50 [cm] x 50 [cm] )の箱のなかに、直方体の鉄塊が ( 20 [cm] x 20 [cm] x 20 [cm] ) が安置してあり、周囲を水が満たしている.

  • 初期温度、 20 [℃] の水が鉄塊に触れて、自然対流し、鉄塊が冷却される.

  • 流体は、すべり境界条件とする.

物性条件

Materials Settings

Target

Parameters

Value

Unit

Description

空気

Heat Conductivity

0.0257

W/m.K

熱伝導度

Heat Capacity

1.005e+3

J/kg.K

比熱

Density

1.166e+0

kg/m3

密度

Heat Expansion Coefficient

3.665e-3

1e-6

体積膨張率(1K毎に(0.29e-3)分の体積増加(線形近似))

Reference Temperature

293.15

K

参照温度 (体積膨張率の測定温度:Boussinesq近似)

viscosity

1.512e-5

Pa.s

粘度

Heat Conductivity

0.602

W/mK

熱伝導度

Heat Capacity

4.18e3

J/kg K

比熱

Density

9.97e3

J/kg K

密度

Reference Temperature

293.15

K

参照温度

Heat Expansion Coefficient

0.29e-3

1/K

体積膨張率 (体積膨張率の測定温度:Boussinesq近似)

viscosity

1.18e-3

Pa s

粘度

SS400

Heat Conductivity

51.6

W/m.K

熱伝導度

Heat Capacity

473.0

J/kg.K

比熱

Density

7.85e+3

kg/m3

密度

定数

gravity

0,0,-1,9.8

m/s2

重力加速度

Stefan Boltzmann

5.68e-8

W/m2K4

Stefan-Boltzmann係数

メッシュ

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

mesh.py ( heatConduction__box_in_water_XYZ3D )
import numpy as np
import os, sys
import gmsh
import nkGmshRoutines.geometrize__fromTable as gft
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"
    dimtagFile = "dat/boundary.json"
    dimtags    = gft.geometrize__fromTable( inpFile=inpFile )
    dimtags    = ldt.load__dimtags( inpFile=dimtagFile, dimtags=dimtags )

    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.050
    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 )
    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(3)
    gmsh.write( "msh/model.msh" )
    gmsh.finalize()

  • geometry.conf

geometry.conf ( heatConduction__box_in_water_XYZ3D )
# <define> @cbox.thick = 5e-3
# <define> @cbox.wx    = 500e-3
# <define> @cbox.wy    = 500e-3
# <define> @cbox.wz    = 500e-3
# <define> @cbox.xMin  = ( -0.5 ) * @cbox.wx
# <define> @cbox.yMin  = ( -0.5 ) * @cbox.wy
# <define> @cbox.zMin  = 0.0

# <define> @water.wx   = @cbox.wx - 2.0*( @cbox.thick )
# <define> @water.wy   = @cbox.wx - 2.0*( @cbox.thick )
# <define> @water.wz   = @cbox.wx -     ( @cbox.thick )
# <define> @water.xMin = ( -0.5 ) * @water.wx
# <define> @water.yMin = ( -0.5 ) * @water.wy
# <define> @water.zMin = @cbox.thick

# <define> @metal.wx   = 200e-3
# <define> @metal.wy   = 200e-3
# <define> @metal.wz   = 200e-3
# <define> @metal.xMin = ( -0.5 ) * @metal.wx
# <define> @metal.yMin = ( -0.5 ) * @metal.wy
# <define> @metal.zMin = @cbox.thick

# <names> key   geometry_type   centering  xc yc zc wx wy wz 
metal     cube  False @metal.xMin @metal.yMin @metal.zMin @metal.wx @metal.wy @metal.wz
water.00  cube  False @water.xMin @water.yMin @water.zMin @water.wx @water.wy @water.wz
cbox      cube  False @cbox.xMin  @cbox.yMin  @cbox.zMin  @cbox.wx  @cbox.wy  @cbox.wz

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

boundary.json ( heatConduction__box_in_water_XYZ3D )
{
    "container.bdr": [ [2,7],[2,8],[2,9],[2,10],[2,19],[2,20] ],
    "metal.bdr"    : [ [2,1],[2,2],[2,3],[2,4],[2,6] ]
}
  • phys.conf

phys.conf ( heatConduction__box_in_water_XYZ3D )
# <names> key	type		dimtags_keys			physNum
metal		      volu		[metal]				    301
water		      volu		[water]				    302
container	    volu		[container]			  303
metal.bdr     surf    [metal.bdr]       201
container.bdr surf    [container.bdr]   202
  • mesh.conf

mesh.conf ( heatConduction__box_in_water_XYZ3D )
# <names> key	physNum	  meshType    resolution1	resolution2	evaluation
metal	  	    301	  constant    0.040		-		-
water		      302	  constant    0.040		-		-
container	    303	  constant    0.040		-		-
metal.bdr     201   constant    0.0     -   -
container.bdr 202   constant    0.0     -   -
  • 生成したメッシュを次に示す.

../../../_images/mesh4.png

Boussinesq近似の確認

シミュレーション設定

  • 流体境界を (i) 鉄塊境界, (ii) 容器境界 の2つに分ける.

    1. 鉄塊境界に温度のDirichlet境界条件 T= 600 [K] = 330 [℃] を課す.

  • 流体としては、 "Normal-Tangential velocity = True" を指定し、滑り壁条件 ( v_\perp=0 ) を課した.

  • Boussinesq近似の挙動を確認する.

boussinesq.sif ( heatConduction__box_in_water_XYZ3D )

!! =============================================================== !!
!! ===  boussinesq.sif ( heatConduction__box_in_water_XYZ3D )  === !!
!! =============================================================== !!

!! -- Boussinesq Approx. Confirmation -- !!
!!
!! * target fluid is "Air"
!! * heat Eqs. in metals is not solved.
!!
!! ------------------------------------- !!

!! ------------------------------------------------- !!
!! --- [1] basics                                --- !!
!! ------------------------------------------------- !!

include "./msh/model/mesh.names"

Header
  CHECK KEYWORDS    Warn
  Mesh DB           "." "msh/model"
  Include Path      ""
  Results Directory "out/"
End

Simulation
  Max Output Level                         = 3
  Coordinate System                        = string "Cartesian"
  Coordinate Mapping(3)                    = 1 2 3

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

  Steady State Max Iterations              = 5
End


Constants
  Stefan Boltzmann                         = 5.6703e-8    !!  -- [ W / m^2 K^4 ] --  !!
  gravity(4)                               = 0 0 -1 9.82  !!  -- [ m/s^2 ] -- !!
End


!! ------------------------------------------------- !!
!! --- [2] solvers / Equations                   --- !!
!! ------------------------------------------------- !!

Solver 1
  Equation                                 = "HeatEquations"
  Procedure                                = "HeatSolve" "HeatSolver"
  Variable                                 = "Temperature"
  Exec Solver                              = "Always"
  Stabilize                                = True
  Bubbles                                  = False
  Optimize Bandwidth                       = True
  Steady State Convergence Tolerance       = 1.0e-5
  Nonlinear System Convergence Tolerance   = 1.0e-3
  Nonlinear System Max Iterations          = 1
  Nonlinear System Newton After Iterations = 3
  Nonlinear System Newton After Tolerance  = 1.0e-3
  Nonlinear System Relaxation Factor       = 1
  Linear System Solver                     = Iterative
  Linear System Iterative Method           = BiCGStab
  Linear System Max Iterations             = 500
  Linear System Convergence Tolerance      = 1.0e-10
  BiCGstabl polynomial degree              = 2
  Linear System Preconditioning            = ILU0
  Linear System ILUT Tolerance             = 1.0e-3
  Linear System Abort Not Converged        = False
  Linear System Residual Output            = 20
  Linear System Precondition Recompute     = 1
End


Solver 2
  Equation                                 = "Navier-Stokes"
  Procedure                                = "FlowSolve" "FlowSolver"
  Variable                                 = Flow Solution[velocity:3,pressure:1]
  Stabilize                                = True
  Bubbles                                  = True
  Optimize Bandwidth                       = True

  Linear System Solver                     = Iterative
  Linear System Scaling                    = Logical False
  Linear System Direct Method              = UMFPACK
  Linear System Iterative Method           = BiCGStab
  Linear System Convergence Tolerance      = 1.0e-6
  Linear System Max Iterations             = 3000
  Linear System Preconditioning            = ILUT

  Nonlinear System Convergence Tolerance   = Real 1.0e-5
  Nonlinear System Max Iterations          = Integer 300
  Nonlinear System Relaxation Factor       = Real 0.7
  Nonlinear System Newton After Iterations = 15
  Nonlinear System Newton After Tolerance  = 1.0e-3

  Steady State Convergence Tolerance       = 1.0e-5

  stabilize                                = True
  Div Discretization                       = True
End

Solver 3
  Exec Solver                              = after saving
  Equation                                 = "Result output"
  Procedure                                = "ResultOutputSolve" "ResultOutputSolver"
  Output File Name                         = "heat"
  Vtu Format                               = Logical True
  Binary Output                            = Logical True
  Scalar Field 1                           = String Temperature
  Scalar Field 2                           = String pressure
  Vector Field 1                           = String velocity
End

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



!! ------------------------------------------------- !!
!! --- [3] Bodies & Materials                    --- !!
!! ------------------------------------------------- !!

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

Material 1
  Name                                     = "Air"
  Heat Conductivity                        = 0.0257     !! W/m.K
  Heat Capacity                            = 1.005e+3   !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Density                                  = 1.166e+0   !! kg/m3
  Heat Expansion Coefficient               = 3.665e-3   !! 1e-6
  viscosity                                = 1.512e-5   !! Pa.s
End


!! ------------------------------------------------- !!
!! --- [5] initial & boundary Condition          --- !!
!! ------------------------------------------------- !!

Initial Condition 1
  Name                                     = "water.initial"
  Temperature                              = 293.15    !! =  20 degree
End                                      


Boundary Condition 1
  Name                                     = "metal.bdr"
  Target Boundaries(1)                     = $metal.bdr
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
  Temperature                              = 600.0
End  

Boundary Condition 2
  Name                                     = "container.bdr"
  Target Boundaries(1)                     = $container.bdr
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
End

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

確認結果

  • 自然対流が生じている.

../../../_images/boussinesq_anime.gif

水中の鉄塊における自然対流冷却

シミュレーション設定

  • Dirichlet条件の場合、空気の動きは見えても、鉄塊内部の冷却の様子はわからない.

  • 解析対象を鉄塊部分の熱伝導も含めて解く.

  • 空気→水へ変更する.

heat.sif ( heatConduction__box_in_water_XYZ3D )

!! ========================================================= !!
!! ===  heat.sif ( heatConduction__box_in_water_XYZ3D )  === !!
!! ========================================================= !!

!! -- Boussinesq Approx. Confirmation -- !!

!! * target fluid is "Air"
!! * heat Eqs. in metals is not solved.

!! ------------------------------------------------- !!
!! --- [1] basics                                --- !!
!! ------------------------------------------------- !!

include "./msh/model/mesh.names"

Header
  CHECK KEYWORDS    Warn
  Mesh DB           "." "msh/model"
  Include Path      ""
  Results Directory "out/"
End

Simulation
  Max Output Level                         = 3
  Coordinate System                        = string "Cartesian"
  Coordinate Mapping(3)                    = 1 2 3

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

  Steady State Max Iterations              = 5
End


Constants
  Stefan Boltzmann                         = 5.6703e-8    !!  -- [ W / m^2 K^4 ] --  !!
  gravity(4)                               = 0 0 -1 9.82  !!  -- [ m/s^2 ] -- !!
End


!! ------------------------------------------------- !!
!! --- [2] solvers / Equations                   --- !!
!! ------------------------------------------------- !!

Solver 1
  Equation                                 = "HeatEquations"
  Procedure                                = "HeatSolve" "HeatSolver"
  Variable                                 = "Temperature"
  Exec Solver                              = "Always"
  Stabilize                                = True
  Bubbles                                  = False
  Optimize Bandwidth                       = True
  Steady State Convergence Tolerance       = 1.0e-5
  Nonlinear System Convergence Tolerance   = 1.0e-3
  Nonlinear System Max Iterations          = 1
  Nonlinear System Newton After Iterations = 3
  Nonlinear System Newton After Tolerance  = 1.0e-3
  Nonlinear System Relaxation Factor       = 1
  Linear System Solver                     = Iterative
  Linear System Iterative Method           = BiCGStab
  Linear System Max Iterations             = 500
  Linear System Convergence Tolerance      = 1.0e-10
  BiCGstabl polynomial degree              = 2
  Linear System Preconditioning            = ILU0
  Linear System ILUT Tolerance             = 1.0e-3
  Linear System Abort Not Converged        = False
  Linear System Residual Output            = 20
  Linear System Precondition Recompute     = 1
End


Solver 2
  Equation                                 = "Navier-Stokes"
  Procedure                                = "FlowSolve" "FlowSolver"
  Variable                                 = Flow Solution[velocity:3,pressure:1]

  Linear System Solver                     = Iterative
  Linear System Scaling                    = Logical False
  Linear System Direct Method              = UMFPACK
  Linear System Iterative Method           = BiCGStab
  Linear System Convergence Tolerance      = 1.0e-6
  Linear System Max Iterations             = 3000
  Linear System Preconditioning            = ILUT

  Nonlinear System Convergence Tolerance   = Real 1.0e-7
  Nonlinear System Max Iterations          = Integer 500
  Nonlinear System Relaxation Factor       = Real 0.7
  Nonlinear System Newton After Iterations = 15
  Nonlinear System Newton After Tolerance  = 1.0e-3

  Steady State Convergence Tolerance       = 1.0e-5

  stabilize                                = True
  Div Discretization                       = False
End

Solver 3
  Exec Solver                              = after saving
  Equation                                 = "Result output"
  Procedure                                = "ResultOutputSolve" "ResultOutputSolver"
  Output File Name                         = "heat"
  Vtu Format                               = Logical True
  Binary Output                            = Logical True
  Scalar Field 1                           = String Temperature
  Scalar Field 2                           = String pressure
  Vector Field 1                           = String velocity
End


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

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



!! ------------------------------------------------- !!
!! --- [3] Bodies & Materials                    --- !!
!! ------------------------------------------------- !!

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

Body 1
  Name                                     = "water"
  Target Bodies(1)                         = $water
  Equation                                 = 1
  Material                                 = 1
  Initial Condition                        = 1
  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                                     = "H2O"
!!   Heat Conductivity                        = 0.602     !! W/m.K
!!   Heat Capacity                            = 4.18e3    !! J/kg.K
!!   Reference Temperature                    = 293.15    !! K
!!   Density                                  = 9.97e+3   !! kg/m3
!!   Heat Expansion Coefficient               = 0.29e-3   !! 1e-6
!!   viscosity                                = 1.18e-3   !! Pa.s
!! End

Material 1
  Name                                     = "Air"
  Heat Conductivity                        = 0.0257     !! W/m.K
  Heat Capacity                            = 1.005e+3   !! J/kg.K
  Reference Temperature                    = 293.15     !! K
  Density                                  = 1.166e+0   !! kg/m3
  Heat Expansion Coefficient               = 3.665e-3   !! 1e-6
  viscosity                                = 1.512e-5   !! Pa.s
End


!! ------------------------------------------------- !!
!! --- [5] initial & boundary Condition          --- !!
!! ------------------------------------------------- !!

!! Initial Condition 1
!!   Name                                     = "SS400.initial"
!!   temperature                              = 673.15    !! = 400 degree
!! End                                      

Initial Condition 1
  Name                                     = "water.initial"
  Temperature                              = 293.15    !! =  20 degree
End                                      

!! Initial Condition 3
!!   Name                                     = "water.initial"
!!   Temperature                              = 293.15    !! =  20 degree
!! End


Boundary Condition 1
  Name                                     = "metal.bdr"
  Target Boundaries(1)                     = $metal.bdr
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
  Temperature                              = 600.0
End  

Boundary Condition 2
  Name                                     = "container.bdr"
  Target Boundaries(1)                     = $container.bdr
  Normal-Tangential Velocity               = True
  velocity 1                               = 0.0
End

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

シミュレーション結果

  • 結果は以下の通り.