Model

class Model(*args)

GetFEM Modelオブジェクト

モデル変数には,変数,状態データ,モデルの説明が格納されます.これには,全体接線行列,RHS,および制約条件が含まれます.モデルには real モデルと complex モデルの2つの種類があります.

Modelオブジェクトの汎用的なコンストラクタ

  • MD = Model('real') 実数が未知数のモデルを作ります.

  • MD = Model('complex') 複素数が未知数のモデルを構築します.

Neumann_term(varname, region)

region のfem変数 varname のNeumann項に対応するアセンブリ文字列を返します.モデルブリックによって宣言されたアセンブリ文字列から推定されます. regionvarname が定義されているメッシュ上の境界領域のインデックスである必要があります.この関数は,すべての体積ブリックが宣言された後で呼び出すようにしてください.ブリックがアセンブリ文字列の宣言を省略した場合はエラーが発生します.

add_Dirichlet_condition_with_Nitsche_method(mim, varname, Neumannterm, datagamma0, region, theta=None, *args)

概要: ind = Model.add_Dirichlet_condition_with_Nitsche_method(self, MeshIm mim, string varname, string Neumannterm, string datagamma0, int region[, scalar theta][, string dataname])

変数 varname とメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります.Neumann項は,高水準汎用アセンブリ言語の表現である(Green公式により得られる)Neumann項の式です.すべての体積ブリックをモデルに追加すると,この項は Model.Neumann_term(varname, region) によって得られます .Dirichlet条件はNitsche法により処理されます. datag はDirichlet条件のオプションのRHSです. datagamma0 はNitsche法のパラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的が強制される標準的な対称法です. theta = -1 は無条件に強制的なskew対称法に対応します. theta = 0 (デフォルト) は非線形問題に対してもNeumann項の二次導関数を必要としない最も単純な方法です.モデル内のブリックインデックスを返します.

add_Dirichlet_condition_with_multipliers(mim, varname, mult_description, region, dataname=None)

変数 varname とメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります.Dirichlet条件は mult_description によって記述される乗数変数で規定されます. mult_description が文字列の場合は,乗数(これは,最初にモデルのメッシュ領域で乗数変数として宣言する必要があります)に対応する変数名と見なされます.有限要素法(mesh_femオブジェクト)の場合は,乗数変数がモデルに追加され,この有限要素法(メッシュ領域 region で制限され,最終的に他の乗数変数との自由度の競合が抑制されます)に基づいて構築されます.整数の場合は,乗数変数がモデルに追加され,その整数の次数の従来の有限要素に基づいて構築されます. dataname はDirichlet条件のオプションのRHSです.定数の場合もあれば,femで記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって異なります.モデル内のブリックインデックスを返します.

add_Dirichlet_condition_with_penalization(mim, varname, coeff, region, dataname=None, mf_mult=None)

変数 varname とメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります.Dirichlet条件はペナルティとともに処理されます.ペナルティ係数は初期値は coeff であり,モデルのデータに追加されます.dataname はDirichlet条件のオプションの右辺です.定数の場合もあれば,有限要素法で記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって異なります. mf_mult はオプションのパラメータで,乗数空間を指定してDirichlet条件を弱めることができます.モデル内のブリックインデックスを返します.

add_Dirichlet_condition_with_simplification(varname, region, dataname=None)

変数 varname とメッシュ領域 region に(単純な)Dirichlet条件を追加します.Dirichlet条件は,(0は対角線の外側,1はマトリックスの対角線,期待値は右側)上の変数の自由度に対応する行を修正することからなる最終的な線形系(非線形問題の接線系)の簡単な後処理によって規定されます.線形システムの対称性は,他のすべてのブリックが対称である場合に維持されます.このブリックは,単純なDirichlet条件(一致境界で宣言された自由度のみが指定されています)用に予約されています.このブリックの還元された自由度への適用は問題があるかもしれません.固有ベクトル有限要素法はサポートされていません. dataname はディリクレ条件のオプションの右側です.これは,定数(この場合は,Lagrange femにのみ適用できます)または varname と同じ有限要素法で記述された(重要な)にすることができます.モデル内のブリックのインデックスを返します.

add_Fourier_Robin_brick(mim, varname, dataexpr, region)

変数 varname に相対的にFourier-Robin項を加えてください.これは以下の形式の形式の弱項に対応します. \(\int (qu).v\)dataexpr はFourier-Robin条件のパラメータ \(q\) です.高水準汎用アセンブリ言語(ただし,モデルのデータである必要がある複雑なバージョンは除きます.)の任意の有効な式を指定できます. region は,項が追加されるメッシュ領域です.モデル内のブリックインデックスを返します.

add_HHO_reconstructed_gradient(transname)

このモデルにHHO法に対する勾配の再構成に対応する基本変換を加えました.名前は基本変換に付けられた名前です.

add_HHO_reconstructed_symmetrized_gradient(transname)

このモデルにHHO法に対する勾配の再構成に対応する基本変換を加えました.名前は基本変換に付けられた名前です.

add_HHO_reconstructed_symmetrized_value(transname)

対称化勾配を用いたHHO法に対する変数の再構成に対応する素変換をモデルに加えます.名前は基本変換に付けられた名前です.

add_HHO_reconstructed_value(transname)

このモデルにHHO法に対する変数の再構成に対応する基本変換を加えました.名前は基本変換に付けられた名前です.

add_HHO_stabilization(transname)

このモデルにHHO安定化演算子に対応する素変換を加えます.名前は,基本変換に付けられた名前です.

add_HHO_symmetrized_stabilization(transname)

対称化勾配を用いたHHO安定化演算子に対応する素変換をモデルに加えます.名前は,基本変換に付けられた名前です.

add_Helmholtz_brick(mim, varname, dataexpr, region=None)

変数 varname に対して相対的にHelmholtz項をモデルに追加します. dataexpr は波数です. region はオプションの項が追加されるメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_Houbolt_scheme(varname)

変数 varname の時間離散化にHoubolt法を適用します.変数の2次時間導関数が多くても存在する場合にのみ有効です.

add_Kirchhoff_Love_Neumann_term_brick(mim, varname, dataname_M, dataname_divM, region)

変数 varname とメッシュ領域 region にKirchhoff-LoveモデルのNeumann項ブリックを追加します. dataname_M は曲げモーメントテンソルを表し, dataname_divM はその発散を表します.モデル内のブリックインデックスを返します.

add_Kirchhoff_Love_plate_brick(mim, varname, dataname_D, dataname_nu, region=None)

変数 varname とメッシュ領域 region にbilaplacianブリックを追加します.これは,項 \(\Delta(D \Delta u)\) を表し,ここで \(D(x)\)dataname_D によって決定される曲げ弾性率です.この項は, dataname_nu をポアソン比とする Kirchhoff-Love プレートモデルに従って部分的に積分されています.モデル内のブリックインデックスを返します.

add_Laplacian_brick(mim, varname, region=None)

変数 varname に相対的にLaplacian項をモデルに追加します(実際には負 \(-\text{div}(\nabla u)\) です).これがベクトル値の変数である場合,Laplace項がコンポーネントごとに追加されます. region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_Mindlin_Reissner_plate_brick(mim, mim_reduced, varname_u3, varname_theta, param_E, param_nu, param_epsilon, param_kappa, variant=None, *args)

概要: ind = Model.add_Mindlin_Reissner_plate_brick(self, MeshIm mim, MeshIm mim_reduced, string varname_u3, string varname_theta , string param_E, string param_nu, string param_epsilon, string param_kappa [,int variant [, int region]])

古典的なReissner-Mindlinプレートモデルに対応する項を追加します.ここで, varname_u3 は横方向の変位, varname_theta は中立面に垂直なファイバーの回転, 'param_E' はヤング率, param_nu はポアソン比, param_epsilon はプレートの厚さ, param_kappa はせん断補正係数です.このブリックは高水準汎用アセンブリ言語を使用しているため,パラメータをこの言語の正規表現にすることができます.3つのバリエーションがあります. variant = 0 は非簡約化された定式化に対応し,その場合には積分法 mim のみが使用されます.実際には,この変形は強いロック現象の影響を受けるため,使用できません. variant = 1 は, mim が回転項に使用され, mim_reduced が横せん断項に使用される縮約積分に対応します. variant = 2 (デフォルト)は,横せん断項の回転RT0要素への投影に対応します.現時点では,(三角形要素のロック現象を除去するだけでは不十分であるため)四角形のみに適用されます.また,高次要素を使用する場合,RT0への投影によって近似の次数が減少することにも注意してください.モデル内のブリックのインデックスを返します.

add_Newmark_scheme(varname, beta, gamma)

変数 varname の時間離散化にtheta法を適用します.変数の2次時間導関数が多くても存在する場合にのみ有効です.

add_Nitsche_contact_with_rigid_obstacle_brick(mim, varname, Neumannterm, dataname_obstacle, gamma0name, region, theta=None, *args)

概要: ind = Model.add_Nitsche_contact_with_rigid_obstacle_brick(self, MeshIm mim, string varname, string Neumannterm, string dataname_obstacle, string gamma0name, int region[, scalar theta[, string dataname_friction_coeff[, string dataname_alpha, string dataname_wt]]])

Coulomb摩擦の有無にかかわらず,接触条件を変数 varname とメッシュ境界 region に追加します.接触条件はNitsche法で規定されています.剛体障害物は,データ dataname_obstacle が障害物までの符号付き距離(有限要素法による補間)で記述されるべきである. gamma0name はNitsche法のパラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的に強制される標準的な対称法に相当する.theta=-1`は無条件に強制的なスキュー対称法に対応する. `theta = 0 はNeumann項の2次導関数を必要としない最も単純な方法です.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定にすることも,有限要素法で定義することもできます.注意: このブリックはNeumann項を持つ偏微分項に対応するすべてのブリックの後にモデルに追加する必要があります.さらに,このブリックはNeumann項を宣言するブリックにのみ適用できます.モデル内のブリックのインデックスを返します.

add_Nitsche_fictitious_domain_contact_brick(mim, varname1, varname2, dataname_d1, dataname_d2, gamma0name, theta=None, *args)

概要: ind = Model.add_Nitsche_fictitious_domain_contact_brick(self, MeshIm mim, string varname1, string varname2, string dataname_d1, string dataname_d2, string gamma0name [, scalar theta[, string dataname_friction_coeff[, string dataname_alpha, string dataname_wt1,string dataname_wt2]]])

仮想ドメイン内の2つのボディ間にクーロン摩擦の有無による接触条件を追加します.接触条件は,Nitsche法では変数 varname_u1 が第1およびスレーブボディに対応し,Nitsche法では変数 varname_u2 が第2およびマスターボディに対応します.接触条件を仮想スレーブ境界上で評価しました.最初のボディはレベルセット dataname_d1 によって記述され,2番目のボディはレベルセット dataname_d2 によって記述されます. gamma0name はNitsche法のパラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的に強制される標準的な対称法に相当する. theta = -1 は無条件に強制的なスキュー対称法に対応する. theta = 0 はNeumann項の2次導関数を必要としない最も単純な方法です.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定にすることも,有限要素法で定義することもできます.注意: このブリックは Neumann 項を持つ偏微分項に対応するすべてのブリックの後にモデルに追加する必要があります.さらに,このブリックはNeumann項を宣言するブリックにのみ適用できます.モデル内のブリックのインデックスを返します.

add_Nitsche_large_sliding_contact_brick_raytracing(unbiased_version, dataname_r, release_distance, dataname_fr=None, *args)

概要: ind = Model.add_Nitsche_large_sliding_contact_brick_raytracing(self, bool unbiased_version, string dataname_r, scalar release_distance[, string dataname_fr[, string dataname_alpha[, int version]]])

Nitsche法に基づいて,摩擦ブリックと有限滑り接触をモデルに追加します.このブリックは自己接触,複数の変形可能体間の接触及び剛体障害物との接触に対処できます.上位レベルの汎用アセンブリを使用します.モデルに raytracing_interpolate_transformation オブジェクトを追加します. "unbiased_version" とは,使用するNitsche法のバージョンのことです.(バイアスされていない,またはバイアスされているもの).スレーブ境界ごとに,この境界上の表示変数の関数として材料法則を定義する必要があります.放出距離は慎重に決定すべきです(一般に平均要素サイズの数倍で,ボディの厚さよりも薄いです).最初は,ブリックは接触境界なしで追加されます.接触境界と剛体には,特別な機能が追加されています. version は非対称バージョンでは0(デフォルト値),より対称なバージョンでは1(摩擦がなくても完全に対称的でない)です.

add_Nitsche_midpoint_contact_with_rigid_obstacle_brick(mim, varname, Neumannterm, Neumannterm_wt, dataname_obstacle, gamma0name, region, theta, dataname_friction_coeff, dataname_alpha, dataname_wt)

実験的なBRICK: 中間点法専用!! Coulomb摩擦の有無にかかわらず,接触条件を変数 varname とメッシュ境界 region に追加します.接触条件はNitsche法で規定されています.剛体障害物は,データ dataname_obstacle が障害物までの符号付き距離(有限要素法による補間)で記述されるべきです. gamma0name はNitsche法のメソッドパラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的に強制される標準的な対称法に相当します. theta = -1 は無条件に強制的なスキュー対称法に対応する. theta = 0 はNeumann項の2次導関数を必要としない最も単純な方法です.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定にすることも,有限要素法で定義することもできます.モデル内のブリックのインデックスを返します.

add_assembly_assignment(dataname, expression, region=None, *args)

概要: Model.add_assembly_assignment(self, string dataname, string expression[, int region[, int order[, int before]]])

アセンブリ時に評価され,im_data型であるデータ dataname に割り当てられる式 expr を追加します.これにより,インスタンスは,他のアセンブリで使用するアセンブリ計算のサブ式を保存できます.たとえば,塑性モデルに塑性歪みを保存するために使用できます.orderは,この割り当てを行うアセンブリの次数(ポテンシャル (0),弱形式 (1),接線システム (2),または各次数(-1))を表します.デフォルト値は1です.before=1の場合,他のアセンブリ項の計算の前に割り当てが実行され,データが残りのアセンブリで中間結果として使用できるようになります(接線システム式の微分が実行されず,データとみなされることに注意してください).before=0(デフォルト)の場合,代入はアセンブリ項の後に行われます.

add_basic_contact_brick(varname_u, multname_n, multname_t=None, *args)

概要: ind = Model.add_basic_contact_brick(self, string varname_u, string multname_n[, string multname_t], string dataname_r, Spmat BN[, Spmat BT, string dataname_friction_coeff][, string dataname_gap[, string dataname_alpha[, int augmented_version[, string dataname_gamma, string dataname_wt]]])

摩擦ブリック付きまたは摩擦ブリックなしの接触をモデルに追加します.Uが,片側拘束が適用される自由度のベクトルである場合,行列 BN はこの拘束が \(B_N U \le 0\) によって定義されるようでなければなりません.摩擦条件は,3つのパラメータ multname_tBTdataname_friction_coeff を追加することで考慮できます.この場合,接線変位は \(B_T U\) であり,行列 BT は, BN\(d-1\) を掛けた数だけの行を持たなければなります.ここで, \(d\) はドメインの次元です.ここでも, dataname_friction_coeff は摩擦係数を表すデータです.各接触条件の値を表すスカラーまたはベクトルです.一方的制約は,その次元が BN の行数に等しくなければならない乗算器 multname_n のおかげで規定される.摩擦条件が追加される場合には,その次元が BT の行数に等しくなければならない乗算器 multname_t によって規定されます.拡張パラメータ r は許容値(Getfemユーザマニュアルを参照)の範囲内で選択するべきです. dataname_gap は初期ギャップを表すオプションのパラメータです.単一の値または値のベクトルを指定できます. dataname_alpha は増大パラメータ(Getfemユーザマニュアルを参照)のためのオプションの均質化パラメータである.パラメータ augmented_version は拡大戦略を示します.1は非対称Alart-Curnier拡大ラグランジアン,2は対称的なもの(接触とCoulomb摩擦の結合を除きます),3は拡大乗算器を用いた非対称方法,4は拡大乗算器と De Saxce 射影を用いた非対称方法です.

add_basic_contact_brick_two_deformable_bodies(varname_u1, varname_u2, multname_n, dataname_r, BN1, BN2, dataname_gap=None, *args)

概要: ind = Model.add_basic_contact_brick_two_deformable_bodies(self, string varname_u1, string varname_u2, string multname_n, string dataname_r, Spmat BN1, Spmat BN2[, string dataname_gap[, string dataname_alpha[, int augmented_version]]])

2つの変形可能モデル間に摩擦なし接触条件を追加します

体.U1,U2が片側拘束が適用される自由度のベクトルである場合,行列 BN1 および BN2 は,この条件が $B_{N1} U_1 B_{N2} U_2 + le gap$ によって定義されるようでなければならない.この制約は,次元が BN の行数に等しくなければならない乗数 multname_n によって規定されます.拡張パラメータ r は許容値(Getfemユーザマニュアルを参照)の範囲内で選択されるべきです. dataname_gap は初期ギャップを表すオプションのパラメータです.単一の値または値のベクトルを指定できます. dataname_alpha は,増大パラメータ(Getfemユーザマニュアルを参照)のためのオプションの均質化パラメータである.パラメータ aug_version は拡大法を示します.1は非対称Alart-Curnier拡大Lagrangian,2は対称Lagrangian,3は拡大乗数を伴う非対称法です.

add_bilaplacian_brick(mim, varname, dataname, region=None)

変数 varname とメッシュ領域 region にbilaplacianブリックを追加します.これは, \(\Delta(D \Delta u)\) という項を表しています.ここで, \(D(x)\)dataname によって決定される係数であり,定数であるか,有限要素法で記述されます.対応する弱形式は \(\int D(x)\Delta u(x) \Delta v(x) dx\) の通りです.モデル内のブリックインデックスを返します.

add_constraint_with_multipliers(varname, multname, B, *args)

概要: ind = Model.add_constraint_with_multipliers(self, string varname, string multname, Spmat B, {vec L | string dataname})

以前にモデルに追加した乗数 multname を使い,変数 varname (は固定サイズ変数でなければなりません)に陽な制約を追加します.制約は \(BU=L\) です.ただし, B は矩形の疎行列です. Model.set_private_matrix() と Model.set_private_rhs() メソッドを使用すると,いつでも拘束を変更できます. L の代わりに dataname が指定された場合,ベクトル L は指定された名前のデータとしてモデル内で定義されます.モデル内のブリックインデックスを返します.

add_constraint_with_penalization(varname, coeff, B, *args)

概要: ind = Model.add_constraint_with_penalization(self, string varname, scalar coeff, Spmat B, {vec L | string dataname})

変数 varname に追加の明示的ペナルティ制約を追加します.制約は以下の通りである. :math`BU=L` と B は矩形の疎行列である. B にはペイン行を含めないでください.そうしないと,接線行列全体が単純になります. Model.set_private_matrix() と Model.set_private_rhs() メソッドを使用すると,いつでも拘束を変更できます.Model.change_penalization_coeff() を使用できます. L の代わりに dataname が指定された場合,ベクトル L は指定された名前のデータとしてモデル内で定義されます.モデル内のブリックインデックスを返します.

add_contact_boundary_to_unbiased_Nitsche_large_sliding_contact_brick(indbrick, mim, region, dispname, lambdaname, wname=None)

接触境界をマスターとスレーブの両方である摩擦ブリックを持つ既存のバイアスされていないNitschelargeスライド接触に追加します.

add_contact_with_rigid_obstacle_brick(mim, varname_u, multname_n, multname_t=None, *args)

概要: ind = Model.add_contact_with_rigid_obstacle_brick(self, MeshIm mim, string varname_u, string multname_n[, string multname_t], string dataname_r[, string dataname_friction_coeff], int region, string obstacle[, int augmented_version])

非推奨機能です.代わりに 'add nodal contact with rigid obstacle brick' を使用してください.

add_data(name, size)

固定サイズのデータをモデルに追加します. sizes は整数(スカラーまたはベクトルデータ用)かテンソルデータの次元のベクトルです. name はデータ名です.

add_elastoplasticity_brick(mim, projname, varname, previous_dep_name, datalambda, datamu, datathreshold, datasigma, region=None)

上位レベルの汎用アセンブリを使用しない古い(廃止予定の)ブリック.等方性材料および準静的モデルの場合,微小変形で変数 varname に対して相対的に非線形弾塑性項をモデルに追加します. projname は使用される投影のタイプです. 'VM' または 'Von Mises' で使用できるのは Von Mises 投影だけです. datasigma はマテリアルの制約を表す変数です. previous_dep_name は,前のタイムステップでのディスプレイスメントを表します.さらに, varname が記述される有限要素法は,K次のmesh_femであり, datasigma は,少なくともK-1次のmesh_femでなければなりません. datalambda および datamu は検討した材料のLame係数である. datathreshold は材料の塑性閾値です.最後の3つの変数は,定数であったり,同じ有限要素法で記述されていたりします. region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_element_extrapolation_transformation(transname, source_mesh, elt_corr)

恒等変換を表すが,多項式外挿によって現在の要素以外の別の要素で式を評価できる特別な補間変換を追加します.仮想領域応用における安定化項に使用されます.配列elt_corは2つのエントリからなる配列で,最初の行は変換に関係する要素を含み,2番目の行は外挿が必要なそれぞれの要素を含みます.要素がelt_corにリストされていない場合,評価は現在の要素に対してのみ行われます.

add_elementary_P0_projection(transname)

投影P0要素に対応する基本変換を追加します.名前は,基本変換に付けられた名前です.

add_elementary_rotated_RT0_projection(transname)

2次元エレメントの回転RT0エレメント上の投影に対応する基本変換をモデルに追加します.名前は,基本変換に付けられた名前です.

add_enriched_Mindlin_Reissner_plate_brick(mim, mim_reduced1, mim_reduced2, varname_ua, varname_theta, varname_u3, varname_theta3, param_E, param_nu, param_epsilon, variant=None, *args)

概要: ind = Model.add_enriched_Mindlin_Reissner_plate_brick(self, MeshIm mim, MeshIm mim_reduced1, MeshIm mim_reduced2, string varname_ua, string varname_theta,string varname_u3, string varname_theta3 , string param_E, string param_nu, string param_epsilon [,int variant [, int region]])

enriched Reissner-Mindlinプレートモデルに対応する項を追加します.ここで, varname_ua は膜変位, varname_u3 は横変位, varname_theta は中間平面に垂直なファイバーの回転,'param_E' はYoung率,param_nu はPoisson比, param_epsilon は板厚です.このブリックは高水準汎用アセンブリ言語を使用しているため,パラメータをこの言語の正規表現にすることができます.4つのバリエーションがあります. variant = 0 は非簡約化された定式化に対応し,その場合には積分法 mim のみが使用されます.実際には,この変形は強いロック現象の影響を受けるため,使用できません. variant = 1 は,mim が回転項に使用され, mim_reduced1 が横せん断項に使用され, mim_reduced2 がピンチ項に使用される縮約積分に対応します. variant = 2 (デフォルト)は,横せん断項の回転RT0要素上への投影と,ピンチ項の縮約積分に対応します.現時点では,(三角形要素のロック現象を除去するだけでは不十分であるため)四角形のみに適用されます.また,高次要素を使用する場合,RT0への投影によって近似の次数が減少することにも注意してください. variant = 3 は,横せん断項の回転RT0要素への投影と,ピンチ項のP0要素への投影とに対応します.現時点では,(三角形要素のロック現象を除去するだけでは不十分であるため)四角形のみに適用されます.また,高次要素を使用する場合,RT0への投影によって近似の次数が減少することにも注意してください.モデル内のブリックのインデックスを返します.

add_explicit_matrix(varname1, varname2, B, issymmetric=None, *args)

概要: ind = Model.add_explicit_matrix(self, string varname1, string varname2, Spmat B[, int issymmetric[, int iscoercive]])

変数 varname1varname2 に対して相対的に接線線形システムに追加される陽な行列を表すブリックを追加します.与えられた行列は varname1 の次元と同じ数の行と varname2 の次元と同じ数の列を持つ必要があります.もし2つの変数が異なっていて,かつ issymmetric が1に設定されているならば,行列の転置もまた接線系に加えられます(デフォルトは0です).項が接線系の保磁力に影響を与えない場合(デフォルトは0です), iscoercive を1に設定します.マトリックスは, Model.set_private_matrix() コマンドで変更できます.モデル内のブリックインデックスを返します.

add_explicit_rhs(varname, L)

変数 varname に対して相対的に接線線形システムの右側に追加される明示的な右辺を表すブリックを追加します.与えられた右辺は varname の次元と同じ大きさでなければなりません.右辺は Model.set_private_rhs() コマンドで変更できます. L の代わりに dataname が指定された場合,ベクトル L は指定された名前のデータとしてモデル内で定義されます.モデル内のブリックインデックスを返します.

add_fem_data(name, mf, sizes=None)

MeshFemにリンクされたモデルにデータを追加します. name はデータ名であり, sizes はオプションのパラメータであり, mf に関して補間的な次元の整数またはベクトルです.

add_fem_variable(name, mf)

MeshFemにリンクされたモデルに変数を追加します. name は変数名です.

add_filtered_fem_variable(name, mf, region)

MeshFemにリンクされたモデルに変数を追加します.変数は,領域上の自由度のみが考慮されるという意味でフィルタリングされます. name は変数名です.

add_finite_strain_elasticity_brick(mim, constitutive_law, varname, params, region=None)

変数 varname を基準にしてモデルに非線形弾性項を追加します. lawname は 'SaintVenant Kirchhoff' , 'Mooney Rivlin' , 'Neo Hookean' , 'Ciarlet Geymonat' ,のいずれかです.対応するバージョンを使用するには,'Mooney Rivlin' および'Neo Hookean' 法の前に 'Compressible' または 'Incompressible' という語を付ける必要があります.これらの法則の圧縮性バージョンには,追加の材料係数が1つ必要です.

重要: 変数が2次元メッシュ上で定義されている場合,平面歪み近似が自動的に使用されます. params は構成則のパラメータのベクトルです.長さは法則によります.これは,定数値の短いベクトル,または可変係数の有限要素法で記述されたベクトルフィールドです.'region' は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.このブリックは,高水準の汎用アセンブリを使用します.モデル内のブリックのインデックスを返します.

add_finite_strain_elastoplasticity_brick(mim, lawname, unknowns_type, varnames=None, *args)

概要: ind = Model.add_finite_strain_elastoplasticity_brick(self, MeshIm mim , string lawname, string unknowns_type [, string varnames, ...] [, string params, ...] [, int region = -1])

有限ひずみ弾塑性ブリックをモデルに追加します.現在のところ,サポートされている法則は, lawname によって "Simo_Miehe" と定義されたものだけです.この法則は, 'DISPLACEMENT_AND_PLASTIC_MULTIPLIER' (整数値1)または 'DISPLACEMENT_AND_PLASTIC_MULTIPLIER_AND_PRESSURE' (整数値3)のいずれかに設定された unknowns_type によって定義される未知変数求解の可能性をサポートします. "Simo_Miehe" の法則では,モデル内で変数として定義する必要がある次の名前の集合を varnames と想定しています.

  • 未知変数として定義されなければならない変位です.

  • 可塑乗数も未知と定義されています.

  • オプションで, 'DISPLACEMENT_AND_PLASTIC_MULTIPLIER_AND_PRESSURE' の混合変位-圧力式の圧力変数を unknowns_type で指定します.

  • 前の時間ステップでの塑性歪みを保持する(スカラー)fem_dataまたはim_dataフィールドの名前です.そして

  • 直前の時間ステップでの右Cauchy-Greenテンソルの逆数の非反復成分をすべて保持するfem_dataまたはim_dataフィールドの名前です(平面歪み2次元問題では4要素ベクトル,3次元問題では6要素ベクトルでなければなりません).

"Simo_Miehe" 法はまた,次の3つのパラメータの集合を params としています.

  • 初期体積弾性率Kの式です,

  • 初期せん断弾性率Gの式です,

  • 硬化変数(降伏限界と硬化変数の値は,それぞれ適切な応力テンソルと歪テンソルのFrobeniusノルムと仮定します)の関数として降伏限界を表示するユーザ定義関数の名前です.

通常, region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_finite_strain_incompressibility_brick(mim, varname, multname_pressure, region=None)

有限ひずみ非圧縮性条件を(有限ひずみ弾性の場合) variable に追加します.multname_pressure は圧力を表す変数である.圧力を記述する有限要素法と主変数との間の入力条件を満たす必要があることに注意してください. region は項が追加されるオプションのメッシュ領域です.指定しない場合はメッシュ全体に追加されます.モデル内のブリックインデックスを返します.p が乗数で,u が変位の場合,"非線形非圧縮性ブリック" で,p*(1-Det(Id(meshdim)+Grad_u)) の高水準汎用アセンブリを追加すると等しいです.

add_generalized_Dirichlet_condition_with_Nitsche_method(mim, varname, Neumannterm, gamma0name, region, theta=None)

変数 varname とメッシュ領域 region にDirichlet条件を追加します.このバージョンはベクトルフィールド用です.これは @f$ Hu = r @f$ という条件を規定しており,ここで H は行列フィールドです.注意: 行列Hの固有値はすべて1または0である必要があります.領域は境界である必要があります. Neumannterm は高水準汎用アセンブリ言語の表現として説明されるNeumann項(Green公式により得られる)の表現です.この項は Model.Neumann_term(varname, region) によってmodelに追加されます.すべてのボリュームブリックがモデルに追加されると.ディリクレ条件はNitsche法で処理されています. dataname はDirichlet条件のオプションの右辺です.これは一定であるか,または有限要素法で記述できます. gamma0name はNitsche法のパラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的に強制される標準的な対称法に相当する. theta = -1 は無条件に強制的なskew対称法に対応します. theta = 0 は非線形問題に対してもNeumann項の2次導関数を必要としない最も単純な方法である. Hname は行列フィールド H に対応するデータである.定数行列であるか,スカラー関数で記述されている必要があります.モデル内のブリックのインデックスを返します.(このブリックは完全にはテストされていせん).

add_generalized_Dirichlet_condition_with_multipliers(mim, varname, mult_description, region, dataname, Hname)

変数 varname とメッシュ領域 region にDirichlet条件を追加します.このバージョンはベクトルフィールド用です.それは H が行列フィールドである条件 \(Hu = r\) を規定します.領域は境界である必要があります.Dirichlet条件は mult_description によって記述される乗数変数で規定される. mult_description が文字列の場合は,乗数(これは,最初にモデルのメッシュ領域で乗数変数として宣言する必要があります)に対応する変数名と見なされます.有限要素法(mesh_femオブジェクト)の場合は,乗数変数がモデルに追加され,この有限要素法(メッシュ領域 region に制限され,最終的に他の乗数変数との自由度の競合が抑制されます)に基づいて構築されます.整数の場合は,乗数変数がモデルに追加され,その整数の次数の従来の有限要素に基づいて構築されます. dataname はDirichlet条件の右側です.定数の場合もあれば,有限要素法で記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって異なります. Hname は行列フィールド H に対応するデータです.モデル内のブリックのインデックスを返します.

add_generalized_Dirichlet_condition_with_penalization(mim, varname, coeff, region, dataname, Hname, mf_mult=None)

変数 varname とメッシュ領域 region にDirichlet条件を追加します.このバージョンはベクトルフィールド用です.これは H が行列フィールドである条件 \(Hu = r\) を規定します.領域は境界である必要があります.ディリクレ状態はペナリゼーションとともに処方される.ペナルティ係数は本来 coeff であり,モデルのデータに追加されます. dataname はDirichlet条件の右辺です.定数の場合もあれば,有限要素法で記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって異なります. Hname は行列フィールド H に対応するデータです.定数行列であるか,スカラー関数で記述されている必要があります. mf_mult はオプションのパラメータで,乗数空間を指定してDirichlet条件を弱めることができます.モデル内のブリックインデックスを返します.

add_generic_elliptic_brick(mim, varname, dataname, region=None)

変数 varname を基準にして,一般的な楕円項をモデルに追加します.楕円項の形状は変数とデータの両方に依存する.これは,以下の項に対応します. \(-\text{div}(a\nabla u)\) ここで, \(a\) はデータ, \(u\) は変数です.データはスカラー,行列,または4次テンソルです.変数はベクトル値であってもなくてもかまいません.データがスカラーまたはマトリックスで,変数がベクトル値の場合,項はコンポーネント単位で追加されます.4次のテンソルデータは,ベクトル値変数にのみ使用できます.データは一定であってもよいし,数分で記述してもよいです.もちろん,データが有限要素法(テンソル場)で記述されるテンソルである場合,データは巨大なベクトルになる可能性があります.行列/テンソルの成分は(blasとの互換性から)データベクトルのFortran(列方向)次数で保存されなければなりません.与えられた行列/テンソルの(仮定されている)対称性は検証されません.これがベクトル値の変数である場合,楕円項はコンポーネントごとに追加されます. region は,項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.高水準汎用アセンブリ言語を使用する実数バージョンでは,モデル変数にも依存しますが, dataname (例えば "1" , "sin(X(1))" ,"Norm(u)" など)は高水準汎用アセンブリ言語の任意の正規表現にすることができます.モデル内のブリックインデックスを返します.

add_im_data(name, mimd)

MeshImdにリンクされたモデルにデータセットを追加します. name はデータ名です.

add_im_variable(name, mimd)

MeshImdにリンクされたモデルに変数を追加します. name は変数名です.

add_initialized_data(name, V, sizes=None)

固定サイズデータを初期化してモデルに追加します. sizes はオプションのパラメータで,データのフォーマットを記述する整数またはベクトルの次元です.デフォルトでは,データはベクトルフィールドbとみなされます. name はデータ名で, V はデータの値です.

add_initialized_fem_data(name, mf, V, sizes=None)

MeshFemにリンクされたモデルにデータを追加します. name はデータ名です.データは V で初期化されます.データはスカラーまたはベクトルフィールドです. sizes はオプションのパラメータであり, `mf`に対する補完次元の整数またはベクトルです.

add_integral_contact_between_nonmatching_meshes_brick(mim, varname_u1, varname_u2, multname, dataname_r, dataname_friction_coeff=None, *args)

概要: ind = Model.add_integral_contact_between_nonmatching_meshes_brick(self, MeshIm mim, string varname_u1, string varname_u2, string multname, string dataname_r [, string dataname_friction_coeff], int region1, int region2 [, int option [, string dataname_alpha [, string dataname_wt1 , string dataname_wt2]]])

摩擦条件の有無にかかわらず,一致しないメッシュ間に接触をモデルに追加します.このブリックは,完全な方法で定義された接触を追加します.これは,連続レベルで定義された拡張Lagrangian定式化(Getfemユーザマニュアルを参照)の直接近似です.この利点は,スケーラビリティの向上です.Newton法の反復回数は,メッシュサイズに多少依存しません.条件は, region1 および region2 に対応する境界上の変数 varname_u1 および varname_u2 に適用されます. multname は,摩擦がない場合の接触応力と,摩擦がある場合の接触および摩擦応力を表す有限要素法の変数である必要があります. multname と`varname_u1` と varname_u2 の間のinf-sup条件が必要です.拡張パラメータ dataname_r は許容値の範囲内で選択する必要があります.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定であるか,または varname_u1 と同じメッシュ上の有限要素法で定義できます. option の値は,非対称 Alart-Curnier 拡張Lagrangian法では1,対称的なものでは2,追加の拡張を伴う非対称 Alart-Curnier 法では3,新しい非対称法では4となります.デフォルト値は1です.摩擦がある場合, dataname_alphadataname_wt1 および dataname_wt2 は,摩擦の進展を解くためのオプションパラメータです.

add_integral_contact_with_rigid_obstacle_brick(mim, varname_u, multname, dataname_obstacle, dataname_r, dataname_friction_coeff=None, *args)

概要: ind = Model.add_integral_contact_with_rigid_obstacle_brick(self, MeshIm mim, string varname_u, string multname, string dataname_obstacle, string dataname_r [, string dataname_friction_coeff], int region [, int option [, string dataname_alpha [, string dataname_wt [, string dataname_gamma [, string dataname_vt]]]]])

摩擦条件の有無にかかわらず,剛体障害物のある接触をモデルに追加します.このブリックは,完全な方法で定義された接触を追加します.これは連続レベルで定義された拡張Lagrangian定式化(Getfemユーザマニュアルを参照)の直接近似である.利点は拡張性が高いことです.ニュートン法の反復回数は,メッシュサイズに多少依存しません.接触条件は region に対応する境界上の変数 varname_u に適用されます.剛体障害物は,データ dataname_obstacle が障害物までの符号付き距離(有限要素法による補間)で記述されるべきである. multname は接触応力を表す有限要素法変数でなければなりません. multnamevarname_u の間のinf-sup条件が必要です.拡張パラメータ dataname_r は,許容値の範囲内で選択する必要があります.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定にすることも,有限要素法で定義することもできます. option の可能な値は,非対称Alart-Curnier拡張Lagrangian法では1,対称的なものでは2,追加の拡張を伴う非対称Alart-Curnier法では3,新しい非対称法では4です.デフォルト値は1です.摩擦がある場合, dataname_alphadataname_wt は摩擦の進化を解決するためのオプションパラメータです.dataname_gamma および dataname_vt は,パラメータに依存するすべり速度を摩擦条件に追加するためのオプションデータを表します.

add_integral_large_sliding_contact_brick_raytracing(dataname_r, release_distance, dataname_fr=None, *args)

概要: ind = Model.add_integral_large_sliding_contact_brick_raytracing(self, string dataname_r, scalar release_distance, [, string dataname_fr[, string dataname_alpha[, int version]]])

摩擦ブリックとの有限すべり接触をモデルに追加します.このブリックは自己接触,複数の変形可能体間の接触及び剛体障害物との接触に対処できます.高水準汎用アセンブリを使用します.モデルに raytracing_interpolate_transformation オブジェクトを追加します.スレーブ境界ごとに乗数変数を定義する必要があります.リリース距離は慎重に決定すべきです(一般に平均要素サイズの数倍で,物体の厚さよりも薄いです).最初は,ブリックは接触境界なしで追加されます.接触境界と剛体には,特別な機能が追加されています. version は,非対称バージョンでは0(デフォルト値),より対称なバージョンでは1(摩擦がなくても完全に対称的でない)です.

add_internal_im_variable(name, mimd)

モデルに変数を追加します.この変数はMeshImdにリンクされており,接線行列のアセンブリ中に縮約されます. name は変数名です.

add_interpolate_transformation_from_expression(transname, source_mesh, target_mesh, expr)

メッシュ source_mesh から式 expr で指定されるメッシュ target_mesh への変換をモデルに追加します.この式は,モデルの変数を含む高水準汎用アセンブリ式に対応します.注意:接線システムの計算では,使用した変数による変換の導関数が考慮されます.ただし,2次微分は実装されていないため,ポテンシャルの定義ではこのような変換は許されていません.

add_isotropic_linearized_elasticity_brick(mim, varname, dataname_lambda, dataname_mu, region=None)

等方性線形化弾性項を,変数 varname に対して相対的にモデルに追加します.dataname_lambdadataname_mu には,Lame係数を含める必要があります. region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_isotropic_linearized_elasticity_pstrain_brick(mim, varname, data_E, data_nu, region=None)

等方性線形化弾性項を,変数 varname に対してモデルに追加します. data_Edata_nu には,それぞれYoung率とPoisson比を含める必要があります. region はオプションの項が追加されるメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.2次元メッシュでは,この項は単純な歪み近似と相関します.3次元メッシュでは,標準モデルに対応します.モデル内のブリックのインデックスを返します.

add_isotropic_linearized_elasticity_pstress_brick(mim, varname, data_E, data_nu, region=None)

等方性線形化弾性項を,変数 varname に対してモデルに追加します.data_Edata_nu には,それぞれYoung率とPoisson比を含める必要があります. region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.2次元メッシュでは,この項は単純応力近似に対応します.3次元メッシュでは,標準モデルに対応します.モデル内のブリックインデックスを返します.

add_linear_generic_assembly_brick(mim, expression, region=None, *args)

概要: ind = Model.add_linear_generic_assembly_brick(self, MeshIm mim, string expression[, int region[, int is_symmetric[, int is_coercive]]])

非推奨.代わりに Model.add_linear_term() を使用します.

add_linear_incompressibility_brick(mim, varname, multname_pressure, region=None, *args)

概要: ind = Model.add_linear_incompressibility_brick(self, MeshIm mim, string varname, string multname_pressure[, int region[, string dataexpr_coeff]])

variable に線形非圧縮性条件を追加します. multname_pressure は圧力を表す変数です.圧力を記述する有限要素法と主変数との間の入力条件を満たす必要があることに注意してください. region は項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます. dataexpr_coeff は,たとえば,非圧縮性に近い弾性に対するオプションのペナルティ係数です.この場合,それはLame係数の逆数である \(\lambda\) です.モデル内のブリックインデックスを返します.

add_linear_term(mim, expression, region=None, *args)

概要: ind = Model.add_linear_term(self, MeshIm mim, string expression[, int region[, int is_symmetric[, int is_coercive]]])

アセンブリ文字列 expr で指定された行列項を追加します.この行列は,領域 region で積分法 mim を使用して組み立てられます.線形であると仮定すると,行列項のみが考慮されます.項を非線形ではなく線形と宣言する利点は,項が1回だけ組み立てられ,残差にアセンブリが必要ないことです.式にいくつかの変数が含まれ,式がポテンシャルまたは一次(すなわち,弱形式の微分ではなく,弱形式を記述する)場合,式はすべての変数に関して微分することに注意してください.項が対称か,強制でないかどうかを指定できます.よくわからない場合は,対称ではなく,強制的ではないと宣言した方がよいでしょう.しかし,一部のソルバー(例えば共役勾配)は非強制問題には使用できません. brickname はブリックのオプション名です.

add_linear_twodomain_term(mim, expression, region, secondary_domain, is_symmetric=None, *args)

概要: ind = Model.add_linear_twodomain_term(self, MeshIm mim, string expression, int region, string secondary_domain[, int is_symmetric[, int is_coercive]])

Modelのような弱形式言語式で与えられる線形項を追加します.add_linear_term()ではなく,二つのドメインの直接の積上での統合の場合には,最初に``mim''と``region''で指定されたものと``secondary_domain''で指定されたものを最初にモデルに宣言しなければなりません.

add_lumped_mass_for_first_order_brick(mim, varname, dataexpr_rho=None, *args)

概要: ind = Model.add_lumped_mass_for_first_order_brick(self, MeshIm mim, string varname[, string dataexpr_rho[, int region]])

変数 varname に対して1次の集中質量項をモデルに追加します.指定された場合,データ dataexpr_rho は密度(省略した場合は1)です. region は,項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックのインデックスを返します.

add_macro(name, expr)

汎用アセンブリ言語の新しいマクロを定義します.名前にはパラメータが含まれます.例えば, name='sp(a,b)' , expr='a.b' は有効な定義です.パラメータのないマクロも定義できます.たとえば, name='x1', expr='X[1]' は有効です. name='grad(u)', expr='Grad_u' という形式も使用できますが,この場合,パラメータ 'u' はマクロの使用時にのみ変数名として許可されます.マクロは,キーワード 'Def' を使用してアセンブリ文字列内に直接定義できます.

add_mass_brick(mim, varname, dataexpr_rho=None, *args)

概要: ind = Model.add_mass_brick(self, MeshIm mim, string varname[, string dataexpr_rho[, int region]])

変数 varname に対して質量項をモデルに追加します.指定された場合,データ dataexpr_rho は密度(省略した場合は1)です. region は,項が追加されるオプションのメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックのインデックスを返します.

add_master_contact_boundary_to_biased_Nitsche_large_sliding_contact_brick(indbrick, mim, region, dispname, wname=None)

摩擦ブリックと既存のバイアスされたNitscheの有限スライド接触にマスター接触境界を追加します.

add_master_contact_boundary_to_large_sliding_contact_brick(indbrick, mim, region, dispname, wname=None)

摩擦ブリックを持つ既存の有限すべり接触にマスター接触境界を追加します.

add_master_contact_boundary_to_projection_transformation(transname, m, dispname, region)

特定の境界 region 上に対応する変位変数 dispname を持つマスター接触境界を,`transname`という既存の投影補間変換に追加します.

add_master_contact_boundary_to_raytracing_transformation(transname, m, dispname, region)

特定の境界 region 上の対応する変位変数 dispname を持つマスター接触境界を transname という既存のレイトレーシング補間変換に追加します.

add_master_slave_contact_boundary_to_large_sliding_contact_brick(indbrick, mim, region, dispname, lambdaname, wname=None)

接触境界を,マスターでもスレーブでもある摩擦ブリックの既存の有限スライド接触に追加します(自己接触は許容します).

add_multiplier(name, mf, primalname, mim=None, region=None)

有限要素法にリンクされた特定の変数を,原始変数に対する乗数として追加します.自由度は gmm::range_basis 関数でフィルタリングされ,乗数と原始変数をリンクするモデルの項に適用されます.これは,主変数に対する線形独立制約のみを保持するためです.境界乗数用に最適化されています.

add_nodal_contact_between_nonmatching_meshes_brick(mim1, mim2=None, *args)

概要: ind = Model.add_nodal_contact_between_nonmatching_meshes_brick(self, MeshIm mim1[, MeshIm mim2], string varname_u1[, string varname_u2], string multname_n[, string multname_t], string dataname_r[, string dataname_fr], int rg1, int rg2[, int slave1, int slave2, int augmented_version])

1つまたは2つの弾性体の2つの面の間に,摩擦条件の有無にかかわらず接触を追加します.条件は,変数 varname_u1 または変数 varname_u1varname_u2 に適用されます.これは,1つまたは2つの個別の変位フィールドが指定されているかどうかによって異なります.整数 rg1 および rg2 は,互いに接触すると予想される領域を表す.単一変位変数の場合, rg1rg2 の両方で定義された領域は,変数`varname_u1` を参照します.2つの変位変数の場合,rg1varname_u1 を参照し,rg2varname_u2 を参照します. multname_t は,rg1 と`rg2` で定義されている領域のうち, "slaves" とされている領域の自由度をサイズとする固定サイズ変数です.接触等価節点垂直力を表します.multname_t は,qdim-1で乗算された multname_n のサイズに対応するサイズの固定サイズ変数でなければなりません.接触等価節点正接(摩擦の)力を表します.拡張パラメータ r は,許容値(弾性体のYoung率に近い値.Getfemのユーザードキュメントを参照)の範囲内で選択されるべきです.パラメータ fr に格納される摩擦係数は,単一の値または multname_n と同じサイズのベクトルのいずれかです.オプションのパラメータ slave1slave2 は, rg1rg2 で定義された領域が "slaves" とみなされるかどうかを宣言します.デフォルトでは, slave1 はtrue, slave2 はfalseです.つまり, 'rg1' はスレーブ面を含み,マスタ面を 'rg2' とします. slave1 と`slave2` のどちらか一方だけをtrueに設定するのが好ましいです.パラメータ augmented_version は拡大戦略を示しmす: 1は非対称のAlart-Curnier拡張Lagrangian,2は対称のもの(接触とCoulomb摩擦の結合を除きます),3は新しい非対称法です.基本的に,このブリックは行列BNとBTとベクトルgapとalphaを計算し,基本的な接触ブリックを呼び出します.

add_nodal_contact_with_rigid_obstacle_brick(mim, varname_u, multname_n, multname_t=None, *args)

概要: ind = Model.add_nodal_contact_with_rigid_obstacle_brick(self, MeshIm mim, string varname_u, string multname_n[, string multname_t], string dataname_r[, string dataname_friction_coeff], int region, string obstacle[, int augmented_version])

摩擦条件の有無にかかわらず,剛体障害物のある接触をモデルに追加します.条件は, region に対応する境界上の変数 varname_u に適用されます.剛障害物は,障害物への符号付き距離である文字列 obstacle で記述します.この文字列は,座標が2次元では 'x', 'y' ,3次元では'x', 'y', 'z' の式である必要があります.たとえば,障害物が \(z \le 0\) に対応する場合,対応する符号付き距離は単に "z" になります. multname_n はサイズが境界 region の自由度数である固定サイズ変数でなければなりません.これは接触等価節点力を表します.摩擦条件を追加するには, multname_t および dataname_friction_coeff パラメータを追加する必要があります.multname_t は固定サイズの変数でなければならず,そのサイズは境界`region` 上の自由度の数に \(d-1\) を乗じたものである.ここで \(d\) は領域の次元です.これは摩擦相当節点力を表します.拡張パラメータ r は許容値(弾性体のYoung率に近い値.Getfemのユーザードキュメントを参照)の範囲内で選択される必要があります.dataname_friction_coeff は摩擦係数です.各接触節点の摩擦係数を表すスカラー値または値のベクトルです.パラメータ augmented_version は拡大戦略を示します: 1は非対称のAlart-Curnier拡張Lagrangian,2は対称(接触とCoulomb摩擦の結合を除く),3は新しい非対称方法です.基本的に,このブリックは行列BNおよびベクトルgapおよびalphaを計算し,基本接触ブリックを呼び出します.

add_nonlinear_elasticity_brick(mim, varname, constitutive_law, dataname, region=None)

変数 varname (deprecated brick,代わりにadd_finite_strain_elasticityを使用してください)に関する非線形弾性項をモデルに追加します.lawname は構成則であり 'SaintVenant Kirchhoff' , 'Mooney Rivlin' , 'neo Hookean' , 'Ciarlet Geymonat' または 'generalized Blatz Ko' のいずれかになります. 'Mooney Rivlin' および 'neo Hookean' 則の名前の前に 'compressible' または 'incompressible' を付けると,対応するバージョンが使用されます.これらの法則の圧縮性バージョンには,追加の材料係数が1つ必要です.既定では,'Mooney Rivlin' 則の非圧縮性バージョンと'neo Hookean' 則の圧縮性バージョンが考慮されます.一般に 'neo Hookean' は'Mooney Rivlin' 則の特殊なケースであり,一つの係数を少なくする必要があります.重要:変数が2次元メッシュ上で定義されている場合,平面歪み近似が自動的に使用されます. dataname は構成則のパラメータのベクトルです.長さは法則によります.これは,定数値の短いベクトル,または可変係数の有限要素法で記述されたベクトルフィールドです. region はオプションの項が追加されるメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.このブリックは,低水準汎用アセンブリを使用します.モデル内のブリックのインデックスを返します.

add_nonlinear_generic_assembly_brick(mim, expression, region=None, *args)

概要: ind = Model.add_nonlinear_generic_assembly_brick(self, MeshIm mim, string expression[, int region[, int is_symmetric[, int is_coercive]]])

非推奨. Model.add_nonlinear_term() を使用します.

add_nonlinear_incompressibility_brick(mim, varname, multname_pressure, region=None)

非線形非圧縮性条件を(有限歪み弾性の場合) variable に追加します. multname_pressure は圧力を表す変数です.圧力を記述する有限要素法と主変数との間の入力条件を満たす必要があることに注意してください. region はオプションの項が追加されるメッシュ領域です.指定しない場合は,メッシュ全体に追加されます.モデル内のブリックインデックスを返します.

add_nonlinear_term(mim, expression, region=None, *args)

概要: ind = Model.add_nonlinear_term(self, MeshIm mim, string expression[, int region[, int is_symmetric[, int is_coercive]]])

アセンブリ文字列 expr で指定された非線形項を追加します.この非線形項は,領域 region と積分法 mim で組み立てられます.この表現は,ポテンシャルまたは弱形式を表すことができます.2次項(すなわち,2次試験関数Test2を含みます)は使用できません.項が対称か,強制かどうかを指定できます.よくわからない場合は,対称ではなく,強制ではないと宣言した方がよいでしょう.しかし,一部のソルバー(例えば共役勾配)は非強制問題には使用できません. brickname はブリックのオプション名です.

add_nonlinear_twodomain_term(mim, expression, region, secondary_domain, is_symmetric=None, *args)

概要: ind = Model.add_nonlinear_twodomain_term(self, MeshIm mim, string expression, int region, string secondary_domain[, int is_symmetric[, int is_coercive]])

弱形式言語の式で与えられる非線形項を追加します. Model.add_nonlinear_term() のようではなく,二つのドメインの直積上での積分の場合には,最初に mimregion で指定されたものと secondary_domain で指定されたものを最初にモデルに宣言しなければなりません.

add_nonmatching_meshes_contact_brick(mim1, mim2=None, *args)

概要: ind = Model.add_nonmatching_meshes_contact_brick(self, MeshIm mim1[, MeshIm mim2], string varname_u1[, string varname_u2], string multname_n[, string multname_t], string dataname_r[, string dataname_fr], int rg1, int rg2[, int slave1, int slave2, int augmented_version])

非推奨機能です.代わりに 'add nodal contact between nonmatching meshes brick' を使用してください.

add_normal_Dirichlet_condition_with_Nitsche_method(mim, varname, Neumannterm, gamma0name, region, theta=None, *args)

概要: ind = Model.add_normal_Dirichlet_condition_with_Nitsche_method(self, MeshIm mim, string varname, string Neumannterm, string gamma0name, int region[, scalar theta][, string dataname])

ベクトル(またはテンソル)値変数 varname の法線コンポーネントとメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります. Neumannterm は,高水準汎用アセンブリ言語の表現として説明されるNeumann項(Green公式により得られる)の表現である.この項は Model.Neumann_term(varname, region) によって得られます.すべてのボリュームブリックがモデルに追加されると,ディリクレ条件はNitscheの方法で処方される. dataname はDirichlet条件のオプションの右辺です.これは一定であるか,または有限要素法で記述できます. gamma0name はNitsche法パラメータです. theta は正または負のスカラー値です. theta = 1gamma0 が小さい場合に条件的に強制される標準的な対称法に相当する. theta = -1 は無条件に強制的なskew対称法に対応します. theta = 0 は非線形問題に対してもNeumann項の2次導関数を必要としない最も単純な方法である.モデル内のブリックのインデックスを返します.(このブリックは完全にはテストされていません)

add_normal_Dirichlet_condition_with_multipliers(mim, varname, mult_description, region, dataname=None)

ベクトル(またはテンソル)値変数 varname の法線コンポーネントとメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります.Dirichlet条件は mult_description によって記述される乗数変数で規定される. mult_description が文字列の場合は,乗数(これは,最初にモデルのメッシュ領域で乗数変数として宣言する必要があります)に対応する変数名と見なされます.有限要素法(mesh_femオブジェクト)の場合は,乗数変数がモデルに追加され,この有限要素法(メッシュ領域 region に制限され,最終的に他の乗数変数との自由度の競合が抑制されます)に基づいて構築されます.整数の場合は,乗数変数がモデルに追加され,その整数の次数の従来の有限要素に基づいて構築されます. dataname はDirichlet条件のオプションの右辺です.定数の場合もあれば,有限要素法で記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって決まります(変数がベクトル値の場合はスカラー,テンソル値の場合はベクトルです).モデル内のブリックのインデックスを返します.

add_normal_Dirichlet_condition_with_penalization(mim, varname, coeff, region, dataname=None, mf_mult=None)

ベクトル(またはテンソル)値変数 varname の法線コンポーネントとメッシュ領域 region にDirichlet条件を追加します.この領域は境界である必要があります.Dirichlet条件はペナリゼーションとともに処理される.ペナルティ係数は最初は coeff であり,モデルのデータに追加されます. dataname はDirichlet条件のオプションの右辺です.定数の場合もあれば,有限要素法で記述される場合もあります.スカラー値またはベクトル値で,Dirichlet条件が指定されている変数によって決まります(変数がベクトル値の場合はスカラー,テンソル値の場合はベクトル). mf_mult はオプションのパラメータで,乗数空間を指定してDirichlet条件を弱めることができます.モデル内のブリックのインデックスを返します.

add_normal_derivative_Dirichlet_condition_with_multipliers(mim, varname, mult_description, region, dataname=None, R_must_be_derivated=None)

変数 varname の法線導関数とメッシュ領域 region (境界である必要がありますにDirichlet条件を追加します).一般的な形式は以下の通りである:math:int partial_n u(x)v(x) = int r(x)v(x) forall v ここで, \(r(x)\) はDirichlet条件(均一条件の場合は0)の右辺であり, \(v\)mult_description で定義される乗数の空間にある. mult_description が文字列の場合は,乗数(これは,最初にモデルのメッシュ領域で乗数変数として宣言する必要があります)に対応する変数名と見なされます.有限要素法(mesh_femオブジェクト)の場合は,乗数変数がモデルに追加され,この有限要素法(メッシュ領域 region に制限され,最終的に他の乗数変数との自由度の競合が抑制されます)に基づいて構築されます.整数の場合は,乗数変数がモデルに追加され,その整数の次数の従来の有限要素に基づいて構築されます. dataname はDirichlet条件の右辺を表すオプションのパラメータです.もし R_must_be_derivatedtrue に設定されていれば, dataname の正規導関数が考慮されます.モデル内のブリックインデックスを返します.

add_normal_derivative_Dirichlet_condition_with_penalization(mim, varname, coeff, region, dataname=None, R_must_be_derivated=None)

変数 varname の法線導関数とメッシュ領域 region (境界である必要があります)にDirichlet条件を追加します.一般的な形式は以下の通りである \(\int \partial_n u(x)v(x) = \int r(x)v(x) \forall v\) ここで \(r(x)\) はDirichlet条件(均一条件の場合は0)の右辺です.ペナルティ係数は最初は coeff であり,モデルのデータに追加されます.これは, Model.change_penalization_coeff() コマンドで変更できます. dataname はDirichlet条件の右辺を表すオプションのパラメータです.もし R_must_be_derivatedtrue に設定されていれば, dataname の正規導関数が考慮されます.モデル内のブリックインデックスを返します.

add_normal_derivative_source_term_brick(mim, varname, dataname, region)

変数 varname とメッシュ領域 region に通常の微分ソース項ブリック \(F = \int b.\partial_n v\) を追加します.

線形システムの右辺を更新します. datanameb を表し, varnamev を表します.モデル内のブリックインデックスを返します.

add_normal_source_term_brick(mim, varname, dataname, region)

境界 region 上の変数 varname にソース項を追加します.この領域は境界である必要があります.ソース項は,データ dataepxpr によって表され,これは,高水準汎用アセンブリ言語(ただし,モデルの宣言されたデータでなければならない複素数バージョンは除きます)の任意の正規表現であり得る.境界への外向き法線単位ベクトルを持つスカラー積が実行されます.このブリックの主な目的は事前処理として法線を持つスカラー積を実行せずにベクトルデータでNeumann条件を表現することです.モデル内のブリックインデックスを返します.

add_penalized_contact_between_nonmatching_meshes_brick(mim, varname_u1, varname_u2, dataname_r, dataname_coeff=None, *args)

概要: ind = Model.add_penalized_contact_between_nonmatching_meshes_brick(self, MeshIm mim, string varname_u1, string varname_u2, string dataname_r [, string dataname_coeff], int region1, int region2 [, int option [, string dataname_lambda, [, string dataname_alpha [, string dataname_wt1, string dataname_wt2]]]])

モデルに,一致しないメッシュ間の摩擦の有無にかかわらずペナルティ接触条件を追加します.条件は, region1 および region2 に対応する境界上の変数 varname_u1 および varname_u2 に適用されます.ペナルティパラメータ dataname_r は,近似の非貫通条件および摩擦条件を規定するのに十分な大きさに選択する必要がありますが,接線システムの条件を過度に悪化させないようにするには大きすぎません.オプションのパラメータ dataname_friction_coeff は摩擦係数で,一定であるか,または varname_u1 と同じメッシュ上の有限要素法で定義できます. dataname_lambda は,optionが2の場合に使用されるオプションのパラメータです.その場合,ペナルティ項はlambda単位でシフトされます(これは,対応する拡張Lagrangian定式化上でのUzawaアルゴリズムの使用を可能にします).接触摩擦の場合, dataname_alphadataname_wt1 および dataname_wt2 は進展摩擦の問題を解決するためのオプションのパラメータです.

add_penalized_contact_with_rigid_obstacle_brick(mim, varname_u, dataname_obstacle, dataname_r, dataname_coeff=None, *args)

概要: ind = Model.add_penalized_contact_with_rigid_obstacle_brick(self, MeshIm mim, string varname_u, string dataname_obstacle, string dataname_r [, string dataname_coeff], int region [, int option, string dataname_lambda, [, string dataname_alpha [, string dataname_wt]]])

摩擦条件の有無にかかわらず,剛体障害物のペナルティ接触をモデルに追加します.条件は region に対応する境界上の変数 varname_u に適用されます.剛体障害物はデータ dataname_obstacle が障害物までの符号付き距離(有限要素法による補間)で記述されるべきである.ペナルティパラメータ dataname_r は,近似の非貫通条件および摩擦条件を規定するのに十分な大きさに選択する必要がありますが,接線システムの条件を過度に悪化させないようにするには大きすぎません. dataname_lambda はオプションが2の場合に使用されるオプションのパラメータです.その場合,ペナルティ項はlambda単位でシフトされます(これは,対応する拡張Lagrangian定式化上でのUzawaアルゴリズムの使用を可能にします).

add_pointwise_constraints_with_given_multipliers(varname, multname, dataname_pt, dataname_unitv=None, *args)

概要: ind = Model.add_pointwise_constraints_with_given_multipliers(self, string varname, string multname, string dataname_pt[, string dataname_unitv] [, string dataname_val])

与えられた乗数 multname を使用して,変数 varname の各点に制約を追加します.条件は,データ dataname_pt に指定された点の集合で規定され,その寸法は点の数とメッシュの寸法の積です.乗数変数は,点数のサイズを固定サイズ変数にする必要があります.変数がベクトルフィールドを表す場合には,データ dataname_unitv を与えなければなりません.これは,次元のベクトルを表すものであり,点の数と単位ベクトルを格納すべきベクトルフィールドの次元を乗算したものです.この場合,指定拘束は,対応する点における変数と対応する単位ベクトルとのスカラー積である.オプションのデータ dataname_val は異なる点で規定される値のベクトルです.このブリックはNeumann問題における剛体変位を消滅させるために特別に設計されています.モデル内のブリックのインデックスを返します.

add_pointwise_constraints_with_multipliers(varname, dataname_pt, dataname_unitv=None, *args)

概要: ind = Model.add_pointwise_constraints_with_multipliers(self, string varname, string dataname_pt[, string dataname_unitv] [, string dataname_val])

multiplierを使用して,変数 varname に点ワイズ制約を追加します.乗数変数が自動的にモデルに追加されます.条件は,データ dataname_pt に指定された点のセットで規定され,その寸法は点の数とメッシュ次元の積である.変数がベクトルフィールドを表す場合には,データ dataname_unitv を与えなければならない.これは,次元のベクトルを表すものであり,点の数と単位ベクトルを格納すべきベクトルフィールドの次元を乗算したものである.この場合,指定拘束は,対応する点における変数と対応する単位ベクトルとのスカラー積である.オプションのデータ dataname_val は,異なる点で規定される値のベクトルです.このブリックは,Neumann問題における剛体変位を消滅させるために特別に設計されている.モデル内のブリックのインデックスを返します.

add_pointwise_constraints_with_penalization(varname, coeff, dataname_pt, dataname_unitv=None, *args)

概要: ind = Model.add_pointwise_constraints_with_penalization(self, string varname, scalar coeff, string dataname_pt[, string dataname_unitv] [, string dataname_val])

ペナルティにより,変数 varname に点ごとの制約を追加します.ペナルティ係数は最初は penalization_coeff で,モデルのデータに追加されます.条件は,データ dataname_pt に指定された点の集合で規定され,その次元は点の数とメッシュの次元の積です.変数がベクトルフィールドを表す場合には,データ dataname_unitv を与えなければなりません.これは,次元のベクトルを表すものであり,点の数と単位ベクトルを格納すべきベクトルフィールドの次元を乗算したものです.この場合,指定拘束は,対応する点における変数と対応する単位ベクトルとのスカラー積です.オプションのデータ dataname_val は異なる点で規定される値のベクトルです.このブリックはNeumann問題における剛体変位を消滅させるために特別に設計されています.モデル内のブリックのインデックスを返します.

add_projection_transformation(transname, release_distance)

transname という名前の投影補間変換を,汎用アセンブリブリックで使用するモデルに追加します.注意: 現時点では,モデルの計算で変換の導関数は考慮されません.

add_raytracing_transformation(transname, release_distance)

汎用アセンブリブリックで使用されるモデルに transname というレイトレーシング補間変換を追加します.注意: 現時点ではモデルの計算で変換の導関数は考慮されません.

add_rigid_obstacle_to_Nitsche_large_sliding_contact_brick(indbrick, expr, N)

摩擦ブリックとの既存の有限すべり接触に,硬質の障害物を追加します. expr は高水準汎用アセンブリ言語(ここで, x はメッシュの現在の点です)を使用した式であり,障害物までの符号付き距離である必要があります. N はメッシュの次元です.

add_rigid_obstacle_to_large_sliding_contact_brick(indbrick, expr, N)

摩擦ブリックとの既存の有限すべり接触に,硬質の障害物を追加します. expr は高水準汎用アセンブリ言語(ここで, x はメッシュの現在の点です)を使用した式であり,障害物までの符号付き距離である必要があります. N はメッシュの次元です.

add_rigid_obstacle_to_projection_transformation(transname, expr, N)

ジオメトリが高水準汎用アセンブリ式 expr の0レベルセットに対応する剛体障害を transname という既存の投影補間変換に追加します.

add_rigid_obstacle_to_raytracing_transformation(transname, expr, N)

ジオメトリが上位レベルの汎用アセンブリ式 expr の0レベルセットに対応する剛体な障害物を `transname`という既存のレイトレーシング補間変換に追加します.

add_slave_contact_boundary_to_biased_Nitsche_large_sliding_contact_brick(indbrick, mim, region, dispname, lambdaname, wname=None)

摩擦ブリックとの既存のバイアスされたNitscheの有限滑り接触にスレーブ接触境界を追加します.

add_slave_contact_boundary_to_large_sliding_contact_brick(indbrick, mim, region, dispname, lambdaname, wname=None)

摩擦ブリックを持つ既存の有限滑り接触にスレーブ接触境界を追加します.

add_slave_contact_boundary_to_projection_transformation(transname, m, dispname, region)

特定の境界 region 上の対応する変位変数 dispname を持つスレーブ接触境界を `transname`という既存の投影補間変換に追加します.

add_slave_contact_boundary_to_raytracing_transformation(transname, m, dispname, region)

特定の境界 region 上の対応するディスプレイスメント変数 dispname を持つスレーブコンタクト境界を, transname という既存のレイトレーシング補間変換に追加します.

add_small_strain_elastoplasticity_brick(mim, lawname, unknowns_type, varnames=None, *args)

概要: ind = Model.add_small_strain_elastoplasticity_brick(self, MeshIm mim, string lawname, string unknowns_type [, string varnames, ...] [, string params, ...] [, string theta = '1' [, string dt = 'timestep']] [, int region = -1])

モデル M に微小歪み塑性項を追加します.これは,微小ひずみ塑性の主要なGetFEMブリックです. lawname は実行された塑性則の名前であり, unknowns_type は,塑性乗数が問題の未知である場合の離散化,または(リターンマッピング法)が次のイテレーションのために格納されるモデルのデータのみの選択を示します.どちらの場合も,乗数は格納されます.varnames は,塑性則に依存する長さの変数名とデータ名の集合です(少なくとも変位,塑性乗数,塑性歪み).paramsはパラメータの式のリストです(少なくとも弾性係数と降伏応力).これらの式には,モデルのデータ名(または変数名)を使用できますが,高水準アセンブリ言語(例えば '1/2', '2+sin(X[0])', '1+Norm(v)' ...)の有効なスカラー式を使用することもできます.パラメータにオプションとして与えられる最後の二つのパラメータは,塑性歪み積分に使用されるスキーム(一般化台形則)のパラメータである theta と時間ステップ dt です.省略された場合の theta のデフォルト値は1であり,これは1次一貫性のある古典的な後方Euler法に対応します. theta=1/2 は二次一貫性のあるCrank-Nicolson法(台形則)に対応します.1/2から1の間の値は有効な値である必要があります. dt のデフォルト値は 'timestep' であり,これは単に(md.set_time_step(dt) により)モデルで定義された時間ステップです.または,任意の式(データ名,定数値...)を使用できます.時間ステップは,1つの反復から次のイテレーションに変更できます. region はメッシュ領域です.

使用可能な塑性則は次のとおりです.

  • 'Prandtl Reuss' (または 'isotropic perfect plasticity' ).硬化のない等方性弾塑性.変数は変位,塑性乗数,塑性歪みです.変位は変数であり,前の時間ステップでの変位(通常は 'u' と 'Previous_u' )に対応する 'Previous_' が先頭に付いた同じ名前の対応するデータを持つ必要があります.塑性乗数には,2つのバージョン(通常は 'xi' と 'Previous_xi' )が必要です.最初のバージョンは, unknowns_type ` が 'DISPLACEMENT_ONLY' または整数値0の場合はデータとして, `unknowns_type がDISPLACEMENT_AND_PLASTIC_MULTIPLIERまたは整数値1の場合は変数として定義されます.塑性歪みは,mesh_femまたは(好ましくは)( mim に対応する)im_dataに保存されている n x n データテンソルフィールドの形である必要があります.データは,最初のLame係数,次の係数(せん断弾性率),および一軸降伏応力です.典型的な呼び出し方法は Model.add_small_strain_elastoplasticity_brick(mim, 'Prandtl Reuss', 0, 'u', 'xi', 'Previous_Ep', 'lambda', 'mu', 'sigma_y', '1', 'timestep'); です.重要:この法則は3次元表現を実装することに注意してください.2次元で使用する場合,式は2次元に単純に転置されます.平面歪み近似については,以下を参照してください.

  • "平面ひずみPrandtl Reuss" (または "平面ひずみ等方性完全塑性")前述の法則と同じですが,平面ひずみ近似に適用されます.2次元でのみ使用できます.

  • "Prandtl Reuss 線形硬化" (または "等方性塑性線硬化").線形等方性および運動硬化を伴う等方性弾塑性. "Prandtl Reuss" 法と比較される追加の変数,累積塑性歪み.塑性歪みと同様に,時間ステップの終了時にのみ保存されるため,単純なデータが必要です(好ましくは, im_data).2つの追加パラメータ: 運動学的硬化係数および等方性.3次元表現のみ.典型的な呼び出しは Model.add_small_strain_elastoplasticity_brick(mim, 'Prandtl Reuss linear hardening', 0, 'u', 'xi', 'Previous_Ep', 'Previous_alpha', 'lambda', 'mu', 'sigma_y', 'H_k', H_i', '1', 'timestep'); です.

  • "平面ひずみPrandtl Reuss線形硬化" (または "平面ひずみ等方性塑性線形硬化").前の法則と同じですが,平面ひずみ近似に適合しています.2次元でのみ使用できます.

塑性流動の離散化および実装された塑性則の詳細についてはGetFEMユーザーマニュアルを参照してください.時間の積分法に関するGetFEMユーザーマニュアル(過渡問題の積分)も参照してください.

重要: small_strain_elastoplasticity_next_iter は各時間ステップの終了時,次の時間ステップの前(と後処理の前: 塑性歪みおよび塑性乗数の値を設定します)に呼び出す必要があります.

add_source_term(mim, expression, region=None)

アセンブリ文字列 expr で指定されたソース項を追加します.これは,領域 region で積分法 mim を使用してアセンブルされます.残差項のみが考慮されます.式にいくつかの変数が含まれていて,式がポテンシャルである場合,式はすべての変数に関して微分されることに注意してください. brickname はブリックのオプション名です.

add_source_term_brick(mim, varname, dataexpr, region=None, *args)

概要: ind = Model.add_source_term_brick(self, MeshIm mim, string varname, string dataexpr[, int region[, string directdataname]])

モデルに変数 varname に比例したソース項を追加します.ソース項は dataexpr で表現します.これは高水準汎用アセンブリ言語の任意の正規表現です(ただし複素数モデルの場合は,宣言されたデータでなければりません). region はオプションで項を追加するメッシュ領域です.追加のオプションデータ directdataname を設定することもできます.対応するデータベクトルは,アセンブリなしで右側に直接追加されます.regionが境界である場合,このブリックでNeumann境界条件を設定することができます.

add_source_term_generic_assembly_brick(mim, expression, region=None)

廃止予定です. Model.add_source_term() を代わりに使用してください.

add_standard_secondary_domain(name, mim, region=-1)

2つのドメインの積を統合するために,弱形式言語の表現で使用できるセカンダリドメインをモデルに追加します. name はセカンダリドメインの名前であり, mim はこのドメイン上の積分法であり, region は積分が実行される領域です.

add_theta_method_for_first_order(varname, theta)

変数 varname の時間離散化にtheta法を適用します.変数の最大でも1次の時間導関数が存在する場合にのみ有効です.

add_theta_method_for_second_order(varname, theta)

変数 varname の時間離散化にtheta法を適用します.変数の2次時間導関数が多くても存在する場合にのみ有効です.

add_twodomain_source_term(mim, expression, region, secondary_domain)

Adds a source term given by a weak form language expression like Model.add_source_term() but for an integration on a direct product of two domains, a first specfied by mim and region and a second one by secondary_domain which has to be declared first into the model.

add_variable(name, sizes)

一定サイズのモデルに変数を追加します. sizes は整数(スカラー変数またはベクトル変数の場合)かテンソル変数の次元のベクトルです. name は変数名です.

assembly(option=None)

すべてのブリックの項を考慮した接線システムのアセンブリです. option を指定する場合, 'build_all' , 'build_rhs' , 'build_matrix' , 'build_rhs_with_internal' , 'build_matrix_condensed' , 'build_all_condensed' とする必要があります.デフォルトでは,接線線形システム全体(行列とrhs)が構築されます.この関数は,独自のソルバを使用して求解する場合に便利です.

brick_list()

モデルのブリックのリストを出力に出力します.

brick_term_rhs(ind_brick, ind_term=None, sym=None, ind_iter=None)

特定の非線形ブリックの項の右側部分にアクセスします.最終的な時間ディスパッチャを考慮しません.まず,右辺のアセンブリを行う必要があります. ind_brick はブリックインデックスです. ind_term はブリック内部の項のインデックスです(デフォルト値は0). sym は2つの異なる変数に作用する対称項のためにの右辺の2番目にアクセスすることです(デフォルトは0です). ind_iter は時間ディスパッチャを使用した場合の反復番号です(デフォルトは0です).

change_penalization_coeff(ind_brick, coeff)

penalizationブリックを使用してDirichlet条件のpenalization係数を変更します.ブリックがこの種類でない場合,この関数の動作は定義されていません.

char()

Modelの(ユニークな)文字列表現を出力します.

これを使用して,2つの異なるModelオブジェクト間の比較を実行できます.この機能は完成予定です.

clear()

モデルをクリアします.

clear_assembly_assignment()

追加されたアセンブリの割り当てをすべて削除します

compute_Von_Mises_or_Tresca(varname, lawname, dataname, mf_vm, version=None)

3次元の非線形弾性のフィールドのVon-Mises応力またはTresca応力を mf_vm で計算します. lawname は構成則で, 'SaintVenant Kirchhoff' , 'Mooney Rivlin' , 'neo Hookean' , 'Ciarlet Geymonat' だ. dataname は構成則のパラメータのベクトルです.長さは法則によります.これは,定数値の短いベクトル,または可変係数の有限要素法で記述されたベクトルフィールドです. version は 'Von_Mises' または 'Tresca' のいずれかです( 'Von_Mises' がデフォルトです).

compute_elastoplasticity_Von_Mises_or_Tresca(datasigma, mf_vm, version=None)

塑性の場のVon-Mises応力またはTresca応力を mf_vm で計算し,ベクトルVに戻します. datasigma はメッシュによってサポートされる応力制約値を含むベクトルです. version は 'Von_Mises' または 'Tresca' ('Von_Mises' がデフォルト)でなければなりません.

compute_finite_strain_elasticity_Von_Mises(lawname, varname, params, mf_vm, region=None)

3次元の非線形弾性のフィールド varname のVon-Mises応力を mf_vm で計算します.lawname は構成則であり,有効な名前である必要があります. params はパラメータ則です.定数値の短いベクトルにすることも,モデルのデータや変数に依存することもできます.高水準汎用アセンブリを使用します.

compute_finite_strain_elastoplasticity_Von_Mises(mim, mf_vm, lawname, unknowns_type, varnames=None, *args)

概要: V = Model.compute_finite_strain_elastoplasticity_Von_Mises(self, MeshIm mim, MeshFem mf_vm, string lawname, string unknowns_type, [, string varnames, ...] [, string params, ...] [, int region = -1])

塑性場のVon-MisesまたはTresca応力を mf_vm で計算し,ベクトルVに返します.最初の入力パラメータは関数 'finite strain elastoplasticity next iter' と同様です.

compute_isotropic_linearized_Von_Mises_or_Tresca(varname, dataname_lambda, dataname_mu, mf_vm, version=None)

Von-Mises 応力またはフィールドのTresca応力(3次元の等方性線形化弾性にのみ有効)を計算します. version は 'Von_Mises' または 'Tresca' ( 'Von_Mises' がデフォルト)でなければなりません.Lame係数によってパラメータ化されます.

compute_isotropic_linearized_Von_Mises_pstrain(varname, data_E, data_nu, mf_vm)

平面ひずみを仮定した3次元または2次元における等方性線形化弾性の変位場のVon-Mises応力を計算します.Young率とPoisson比によってパラメータ化されます.

compute_isotropic_linearized_Von_Mises_pstress(varname, data_E, data_nu, mf_vm)

平面応力を仮定した3次元または2次元における等方性線形化弾性の変位場のVon-Mises応力を計算します.Young率とPoisson比によってパラメータ化されます.

compute_plastic_part(mim, mf_pl, varname, previous_dep_name, projname, datalambda, datamu, datathreshold, datasigma)

プラスチック成形品を mf_pl で計算し,ベクトルVに返します. datasigma はメッシュでサポートされる応力拘束値を含むベクトルです.

compute_second_Piola_Kirchhoff_tensor(varname, lawname, dataname, mf_sigma)

3次元非線形弾性の場の第2Piola Kirchhoff応力テンソルのmf_sigmaを計算します. lawname は構成則で 'SaintVenant Kirchhoff' , 'Mooney Rivlin' , 'neo Hookean' ,あるいは 'Ciarlet Geymonat'となります. dataname は構成則のパラメータのベクトルです.長さは法則によります.これは,定数値の短いベクトル,または可変係数の有限要素法で記述されたベクトルフィールドです.

contact_brick_set_BN(indbrick, BN)

基本的な接触/摩擦ブリックのBNマトリックスの設定に使用できます.

contact_brick_set_BT(indbrick, BT)

摩擦ブリックとの基本接触のBTマトリックスの設定に使用できます.

define_variable_group(name, varname=None, *args)

概要: Model.define_variable_group(self, string name[, string varname, ...])

補間 (主にレイトレーシング補間変換の変数のグループを定義します.

del_macro(name)

汎用アセンブリ言語用に定義済みのマクロを削除します.

delete_brick(ind_brick)

モデルから変数またはデータを削除します.

delete_variable(name)

モデルから変数またはデータを削除します.

disable_bricks(bricks_indices)

ブリックを無効にします(ブリックは接線線形システムの構築に関与しなくなります).

disable_variable(varname)

求解のため変数(とその添付乗数)を使用不可にします.次の求解は,残りの変数に対してのみ実行されます.これにより,モデルの異なる部分を別々に求解できます.変数の強連成がある場合は,固定小数点方式を使用できます.

displacement_group_name_of_Nitsche_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

displacement_group_name_of_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

display()

Modelオブジェクトの概要が表示されます.

elastoplasticity_next_iter(mim, varname, previous_dep_name, projname, datalambda, datamu, datathreshold, datasigma)

古い(廃止された)弾塑性ブリックとともに使用して,イテレーションから次のイテレーションに渡します.次のイテレーションの応力拘束sigmaを計算して保存します. 'mim' は,計算に使用する積分法です. 'varname' がこの問題の主要な変数です. 'previous_dep_name' は前の時間ステップでの変位を表します. 'projname' は使用する投影のタイプです.現時点では 'Von Mises' または 'VM' でなければなりません. 'datalambda' および 'datamu' は,材料のLame係数です. 'datasigma' は新しい応力拘束値を格納するベクトルです.

enable_bricks(bricks_indices)

無効なブリックを有効にします.

enable_variable(varname)

使用不可の変数を使用可能にします(付属の乗数).

finite_strain_elastoplasticity_next_iter(mim, lawname, unknowns_type, varnames=None, *args)

概要: Model.finite_strain_elastoplasticity_next_iter(self, MeshIm mim, string lawname, string unknowns_type, [, string varnames, ...] [, string params, ...] [, int region = -1])

有限ひずみ塑性ブリックの時間ステップから別のステップへの移行を可能にする関数です.パラメータは, add_finite_strain_elastoplasticity_brick とまったく同じでなければなりません.説明については,この関数のドキュメントを参照してください.基本的に,このブリックは塑性歪みと塑性乗数を計算し,次のステップのために保存します.現在実装されている唯一のSimo-Mieheの法則に対して,この関数は varnames の最後の2つのエントリで定義されている状態変数を更新し, varnames の2番目のエントリとして与えられる塑性乗数フィールドをリセットします.

first_iter()

時間積分スキームの最初の反復の前に実行します.

from_variables()

モデルの変数を連結したモデルのすべての自由度(独自のソルバを使用して求解するのに便利)のベクトルを返します.

get_time()

現在時刻に対応するデータ t の値を指定します.

get_time_step()

タイムステップの値を指定します.

interpolation(expr, *args)

概要: V = Model.interpolation(self, string expr, {MeshFem mf | MeshImd mimd | vec pts, Mesh m}[, int region[, int extrapolation[, int rg_source]]])

mesh_fem mf または mesh_im_data mimd またはメッシュ m 上の点の集合 pts を基準にして特定の式を補間します.式は,モデルの変数とデータへの参照を含む可能性があります,高水準汎用アセンブリ言語に従って有効である必要があります.

オプション 外挿rg_source は,点 pts の集合に関する補間のためのものです.

interval_of_variable(varname)

モデルの線形システムの変数 varname の間隔を指定します.

is_complex()

モデルが実数の場合は0を返し,複素数の場合は1を返します.

list_residuals()

モデルに含まれるすべての項に対応する残差を出力します.

local_projection(mim, expr, mf, region=None)

mesh_fem mf を基準にして,式の要素ごとのL2投影を作成します.このmesh_femは不連続である必要があります.式は,モデルの変数とデータへの参照を含む可能性があります,上位レベルの汎用アセンブリ言語に従い有効である必要があります.

matrix_term(ind_brick, ind_term)

存在する場合,ブリック ind_brick のマトリックス項 ind_term を返します.

memsize()

モデルが使用するメモリ量(バイト単位)の概算値を返します.

mesh_fem_of_variable(name)

変数またはデータの mesh_fem へのアクセスを許可します.

mult_varname_Dirichlet(ind_brick)

Dirichletブリックの乗数変数の名前を指定します.ブリックが乗数ブリックを持つDirichlet条件ではない場合,この関数の動作は未定義です.

nbdof()

モデルの自由度の総数を返します.

next_iter()

時間積分スキームの各反復の終了時に実行されます.

perform_init_time_derivative(ddt)

この関数を呼び出すことにより,時間積分スキームによって必要とされる導関数に対応するデータを初期化するためには,次の求解が(非常に)小さな時間ステップ ddt (主に時間問題では1次の初期時間微分,時間問題では2次の2次時間微分)の解を計算することを示します.次の求解では,変数の値は変更されません.

resize_variable(name, sizes)

モデルの一定サイズ変数のサイズを変更します. sizes は(スカラー変数またはベクトル変数の場合)整数かテンソル変数の次元のベクトルです. name は変数名です.

rhs()

接線問題の右辺を返します.

set_element_extrapolation_correspondence(transname, elt_corr)

要素外挿補間変換の対応マップを変更します.

set_private_matrix(indbrick, B)

内部に疎行列がある特定のブリック(陽なブリック: '拘束ブリック' と '陽な行列ブリック' )では,この行列を設定します.

set_private_rhs(indbrick, B)

内部右辺ベクトル(陽なブリック: '拘束ブリック' と '陽なブリック')を持つ特定のブリックでは,このrhsを設定します.

set_time(t)

現在時刻 t に対応するデータの値を t に設定します.

set_time_step(dt)

時間ステップの値を dt に設定します.この値は,すべてのワンステップスキームでステップから別のステップに変更できます(すなわち現時点では,提案されているすべての時間積分法).

set_variable(name, V)

変数またはデータの値を設定します. name はデータ名です.

shift_variables_for_time_integration()

モデルの変数を,時間積分スキームの前の時間ステップの他の値に対応するデータにシフトするために使用される関数です.時間積分スキームが宣言されている変数ごとに,スキームが呼び出されてシフトが実行されます.この関数は,2つのタイムステップの間に呼び出す必要があります.

sliding_data_group_name_of_Nitsche_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

sliding_data_group_name_of_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

small_strain_elastoplasticity_Von_Mises(mim, mf_vm, lawname, unknowns_type, varnames=None, *args)

概要: V = Model.small_strain_elastoplasticity_Von_Mises(self, MeshIm mim, MeshFem mf_vm, string lawname, string unknowns_type [, string varnames, ...] [, string params, ...] [, string theta = '1' [, string dt = 'timestep']] [, int region])

この関数は, mf_vm で近似された小さな歪みの弾塑性項に関するフォンミーゼス応力場を計算し,その結果を VM に格納します.その他のパラメータはすべて, add_small_strain_elastoplasticity_brick とまったく同じでなければなりません.この関数を呼び出す前に, small_strain_elastoplasticity_next_iter を呼び出す必要があります.

small_strain_elastoplasticity_next_iter(mim, lawname, unknowns_type, varnames=None, *args)

概要: Model.small_strain_elastoplasticity_next_iter(self, MeshIm mim, string lawname, string unknowns_type [, string varnames, ...] [, string params, ...] [, string theta = '1' [, string dt = 'timestep']] [, int region = -1])

微小ひずみ塑性ブリックの時間ステップ間を通過できる関数.パラメータは,add_small_strain_elastoplasticity_brick とまったく同じでなければなりません.説明については,この関数のドキュメントを参照してください.基本的に,このブリックは塑性歪みと塑性乗数を計算し,次のステップのために保存します.さらに,計算されたディスプレイスメントが,前のタイムステップのディスプレイスメントを格納するデータにコピーされます(通常は 'u' から 'Previous_u' ).これは, compute_small_strain_elastoplasticity_Von_Mises を使用する前に呼び出す必要があります.

solve(*args)

概要: (nbit, converged) = Model.solve(self[, ...])

標準のgetfemソルバーを実行します.

必要に応じて,独自のソルバを使用できるようにする必要があります( Model.tangent_matrix() などで接線行列とその右辺を得ることができます).

さまざまなオプションを指定できます.

  • 'noisy' または 'very_noisy'

    ソルバが進行状況を示す情報が表示します(残差値など).

  • 'max_iter', int NIT

    最大反復回数を設定します.

  • 'max_res', @float RES

    目標残差値を設定します.

  • 'diverged_res', @float RES

    反復法が発散したと見なされる残差のしきい値を設定します(デフォルトは1e200です).

  • 'lsolver', string SOLVER_NAME

    線形システム(デフォルト値は 'auto' で,getfemが自分で選択できるようになっています)に使用するソルバーを明示的に選択します.指定できる値は 'superlu' ,(サポートされている場合) 'mumps' , 'cg/ildlt' , 'gmres/ilut' ,および 'gmres/ilut' です.

  • 'lsearch', string LINE_SEARCH_NAME

    線形システムで使用される線検索方法を明示的に選択します(デフォルト値は 'default' です).指定できる値は, 'simplest' , 'systematic' , 'quadratic' ,または 'basic' です.

    反復法が使用される場合,反復回数を返します.

    変数のサブセット(無効にされた変数はデータとみなされます)に関してのみ問題を求解するために,いくつかの変数( Model.disable_variable() を参照してください)を無効にすることが可能であることに注意してください.例えば,全体Newton法を固定小数点法で置き換えることができます.

tangent_matrix()

モデルに格納されている接線行列を返します.

test_tangent_matrix(EPS=None, *args)

概要: Model.test_tangent_matrix(self[, scalar EPS[, int NB[, scalar scale]]])

ランダムな位置と方向(新しく作成されたブリックをテストするのに便利)で,接線行列の整合性をテストします.EPSは微分の差分計算のための小パラメータの値であり,ランダム方向(デフォルトは1E-6)です. NN はテストの数です(デフォルトは100です). scale は現在位置の周囲のランダムな位置(デフォルトは1で,0が許容値です)のパラメータです.ランダム位置の各自由度は [current-scale, current+scale] の範囲内で選択されます.

test_tangent_matrix_term(varname1, varname2, EPS=None, *args)

概要: Model.test_tangent_matrix_term(self, string varname1, string varname2[, scalar EPS[, int NB[, scalar scale]]])

ランダムな位置と方向(新しく作成されたブリックをテストするのに便利)で,接線行列の一部の整合性をテストします.インクリメントは,変数 varname2 に対してのみ行われ, varname1 に対応する残差の部分でテストされます.これは,接線行列の項 (varname1, varname2) のみがテストされることを意味します.EPSは微分の差分計算のための微小パラメータの値(デフォルトは1E-6)であり,ランダム方向です. NN はテストの数です(デフォルトは100です). scale は現在位置の周囲のランダムな位置(デフォルトは1で,0が許容値です)のパラメータである.ランダム位置の各自由度は [current-scale, current+scale] の範囲内で選択されます.

to_variables(V)

ベクトル V でモデルの変数の値を設定します.通常,ベクトル V は接線線形システム(独自のソルバを使用して問題を解決するのに便利)のソルバの結果です.

transformation_name_of_Nitsche_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

transformation_name_of_large_sliding_contact_brick(indbrick)

既存の有限すべり接触ブリックの滑りデータに対応する変数のグループの名前を指定します.

variable(name)

変数またはデータの値を指定します.

variable_list()

モデルの変数と定数のリストを出力に出力します.