コンテンツにスキップ

ユーザーサブルーチン

ユーザーサブルーチン

ユーザーがFrontISTRの機能をプログラミングにより拡張するためのインターフェースを提供する。 これらのインターフェースは、基本的にサブルーチンヘッダを含むFORTRANサブルーチンで、入出力変数の記述とこれらの変数のための宣言文である。 ルーチンの主要部は、ユーザーによって書かなければならない。

FrontISTRは以下のユーザサブルーチンインターフェースを提供している。

ユーザー定義材料の入力

ユーザー定義材料を使用する場合、最大100のユーザー定義材料定数が使用可能である。 材料定数の入力は以下のように、制御データファイル内の1行10数値、最大10行まで入力可能である。

2行目~最大10行目

v1, v2, v3, v4, v5, v6, v7, v8, v9, v10

弾塑性変形に関わるサブルーチン (uyield.f90)

弾塑性剛性マトリクスおよび応力のreturn mappingを計算するためのサブルーチンを提供している。 ユーザー定義降伏関数を利用する場合、まず入力ファイルに!PLASTIC, YIELD=USERを設定して必要な材料定数を入力し、 次に関数uElastoPlasticNumStatusとサブルーチンuElastoPlasticMatrixおよびuBackwardEuler を作成する必要がある。

(1) 実数の状態変数の数を返す関数

integer function uElastoPlasticNumStatus( matl )
    real(kind=kreal),   intent(in)  :: matl(:)
  • matl: 材料定数を保存する配列(1-100 : システム定義の材料定数、101-200 : ユーザー定義の材料定数)

(2) 弾塑性剛性マトリクスの計算サブルーチン

  subroutine uElastoPlasticMatrix( matl, stress, istat, fstat, plstrain, D, temp, hdflag )
    real(kind=kreal),   intent(in)  :: matl(:)
    real(kind=kreal),   intent(in)  :: stress(6)
    integer(kind=kint), intent(in)  :: istat
    real(kind=kreal),   intent(in)  :: fstat(:)
    real(kind=kreal),   intent(in)  :: plstrain
    real(kind=kreal),   intent(out) :: D(:,:)
    real(kind=kreal),   intent(in)  :: temp
    integer(kind=kint), intent(in)  :: hdflag
  • matl: 材料定数を保存する配列(1-100 : システム定義の材料定数、101-200 : ユーザー定義の材料定数)
  • stress: 2nd Piola-Kirchhoff応力
  • istat: 整数の状態変数
  • fstat: 実数の状態変数の配列
  • plstrain: 現在のサブステップ開始時点の塑性ひずみ
  • D: 弾塑性マトリクス
  • temp: 温度
  • hdflag: 全成分(0)、偏差成分のみ(1)、体積成分のみ(2)を計算

(3) 応力のReturn mapping計算サブルーチン

  subroutine uBackwardEuler( matl, stress, plstrain, istat, fstat, temp, hdflag )
    real(kind=kreal),   intent(in)    :: matl(:)
    real(kind=kreal),   intent(inout) :: stress(6)
    real(kind=kreal),   intent(in)    :: plstrain
    integer(kind=kint), intent(inout) :: istat
    real(kind=kreal),   intent(inout) :: fstat(:)
    real(kind=kreal),   intent(in)    :: temp
    integer(kind=kint), intent(in)    :: hdflag
  • matl: 材料定数を保存する配列(1-100 : システム定義の材料定数、101-200 : ユーザー定義の材料定数)
  • stress: trial stress弾性変形を仮定し得られた2nd Piola-Kirchhoff応力
  • plstrain: 現在のサブステップ開始時点の塑性ひずみ
  • istat: 整数の状態変数
  • fstat: 実数の状態変数の配列
  • temp: 温度
  • hdflag: 全成分(0)、偏差成分のみ(1)、体積成分のみ(2)を計算

弾性変形に関わるサブルーチン (uelastic.f90)

弾性および超弾性問題の弾性剛性マトリクスおよび応力の更新計算をするためのサブルーチンを提供している。ユーザー弾性または超弾性構成式を利用する場合、まず入力ファイルに!ELASTIC, TYPE=USERまたは!HYPERELASTIC, TYPE=USERを設定して必要な材料定数を入力し、次にサブルーチンuElasticMatrixおよびuElasticUpdateを作成する必要がある。

(1) 弾性剛性マトリクスの計算サブルーチン

subroutine uElasticMatrix( matl, strain, D )
    REAL(KIND=kreal), INTENT(IN) :: matl(:)
    REAL(KIND=kreal), INTENT(IN) :: strain(6)
    REAL(KIND=kreal), INTENT(OUT) :: D(6,6)
  • matl: 材料定数を保存する配列(最大100)
  • strain: Green-Lagrangeひずみ
  • D: 弾性マトリクス

(2) 応力の計算サブルーチン

subroutine uElasticUpdate ( matl, strain, stress )
    REAL(KIND=kreal), INTENT(IN) :: matl(:)
    REAL(KIND=kreal), INTENT(IN) :: strain(6)
    REAL(KIND=kreal), INTENT(OUT) :: stress(6)
  • matl: 材料定数を保存する配列(最大100)
  • strain: Green-Lagrangeひずみ
  • stress: 応力

ユーザー定義材料に関わるサブルーチン (umat.f)

弾性、超弾性、弾塑性材に拘らず一般的な材料の変形解析のインターフェースを提供する。 ユーザー定義材料を利用する場合、まず入力ファイルに!USER_MATERIALを設定して必要な材料定数を入力し、次にサブルーチンuMatlMatrixおよびuUpdateを作成する必要がある。

(1) 剛性マトリクスの計算サブルーチン

subroutine uMatlMatrix( mname, matl, strain, stress, fstat, D, dtime, ttime, temperature )
        character(len=*), intent(in)  :: mname
        real(kind=kreal), intent(in)  :: matl(:)
        real(kind=kreal), intent(in)  :: strain(6)
        real(kind=kreal), intent(in)  :: stress(6)
        real(kind=kreal), intent(in)  :: fstat(:)
        real(kind=kreal), intent(out) :: D(:,:)
        real(kind=kreal), intent(in)  :: dtime
        real(kind=kreal), intent(in)  :: ttime
        real(kind=kreal), optional    :: temperature
  • mname: 材料名
  • matl: 材料定数を保存する配列(最大100)
  • strain: Green-Lagrangeひずみ
  • stress: 2nd Piola-Kirchhoff応力
  • fstat: 状態変数
  • D: 構成式
  • dtime: 時間増分
  • ttime: 現在の時間増分開始時点のトータル時刻
  • temperature: 温度

(2) ひずみおよび応力の更新計算サブルーチン

subroutine uUpdate(  mname, matl, strain, stress, fstat, dtime, ttime, temperature )
        character(len=\*), intent(in)    :: mname
        real(kind=kreal), intent(in)    :: matl(:)
        real(kind=kreal), intent(in)    :: strain(6)
        real(kind=kreal), intent(inout) :: stress(6)
        real(kind=kreal), intent(inout) :: fstat(:)
        real(kind=kreal), intent(in)    :: dtime
        real(kind=kreal), intent(in)    :: ttime
        real(kind=kreal), optional      :: temperature
  • mname: 材料名
  • matl: 材料定数を保存する配列(最大100)
  • strain: ひずみ
  • stress: 2nd Piola-Kirchhoff応力
  • fstat: 状態変数
  • dtime: 時間増分
  • ttime: 現在の時間増分開始時点のトータル時刻
  • temperature: 温度

ユーザー定義外部荷重の処理サブルーチン (uload.f)

ユーザー定義外部荷重を処理するインターフェースを提供する。

ユーザー定義外部荷重を利用するため、まず外部荷重を定義するための数値構造tULoadを定義し、入力ファイルの!ULOAD, FILE=ファイル名を利用してファイル名で指定されたファイルからその定義を読み込む。 その後、以下のインターフェースを利用して、外部荷重を組み込む。

(1) 外部荷重の読み込みサブルーチン

integer function ureadload( fname )
    character(len=\*), intent(in) :: fname
  • fname: 外部ファイル名。このファイルからユーザー定義外部荷重を読み込む。

(2) 外部荷重を全体荷重ベクトルへ組み込むサブルーチン

subroutine uloading( cstep, factor, exForce )
    integer, INTENT(IN) :: cstep
    REAL(KIND=kreal), INTENT(IN) :: factor
    REAL(KIND=kreal), INTENT(INOUT) :: exForce(:)
  • cstep: 現時点の解析ステップ数
  • factor: 現ステップの荷重係数
  • exForce: 全体荷重ベクトル

(3) 残差応力の計算サブルーチン

subroutine uResidual( cstep, factor, residual )
    integer, INTENT(IN) :: cstep
    REAL(KIND=kreal), INTENT(IN) :: factor
    REAL(KIND=kreal), INTENT(INOUT) :: residual(:)
  • cstep: 現時点の解析ステップ数
  • factor: 現ステップの荷重係数
  • residual: 全体残差力ベクトル