本ライブラリルーチンで用いた各反復法アルゴリズムについて紹介する. 以下では,nnの実行列Aを係数行列に持つ連立一次方程式 を解くものとする.ここでは,適当な前処理行列を,初期値ベクトルをとする.
共役勾配法(CG法)は対称行列向けのクリロフ部分空間法である.
安定化双共役勾配法(Bi-CGSTAB)は非対称行列向けのクリロフ部分空間法である.
一般化最小残差法(GMRES(m))は非対称行列向けのクリロフ部分空間法である.リスタート周期m回の反復計算で解が収束しない場合はm回の反復で得られた解を初期値として反復を再開する.PARCELでは直交化手法として古典グラムシュミット法と修正グラムシュミット法が整備されている.
はを要素に持つ上三角行列,はベクトルの最初のi成分からなるベクトルである.
通信削減型一般化最小残差法(CA-GMRES(s,t))は非対称行列向けのクリロフ部分空間法である. GMRES(m)のs回反復に相当する計算を1反復で行う. t回の反復で収束しない場合はt回の反復で得られた近似解を初期解として反復を再開する. m = t sがGMRES(m)のリスタート周期に相当する. 数学的にはCA-GMRES(s,t)とGMRES(m)の収束性は同等であり,通信削減型のQR分解を適用することでGMRES(m)よりも集団通信の回数を削減する. ただし,一度にs本の基底ベクトルを生成するため,sが過大な場合には,丸め誤差によって基底ベクトルの直交性が崩れることでGMRESよりも収束性が悪くなる.PARCELでは基底ベクトルの直交性改善のために単基底に加えてNewton基底が整備されている.
ここで,は単位ベクトル, はのs行s列要素, は反復の過程で得られるヘッセンベルグ行列から求めた固有値である.
はを要素に持つ行列.
はを要素に持つ行列.
チェビシェフ基底共役勾配法(CBCG)は対称行列向けのクリロフ部分空間法である.CG法のk反復に相当する計算を1反復で行うことで集団通信の回数を削減する.チェビシェフ基底を用いて基底ベクトルを生成するため行列の最大・最小固有値を必要とする.PARCELでは固有値計算法としてべき乗法と通信削減型アーノルディ法が整備されている.
,はの最大及び最小固有値.
一般に反復解法では,となるように前処理行列を選ぶと,収束に必要な反復回数が減少する. アルゴリズムに見られる 「」型の計算は 「 を について近似的に解く」ことに相当する. この近似の精度を上げれば反復回数が減少するが,逆にこの部分の計算のオーバーヘッドが増加する. 並列計算でない場合には,この近似計算に不完全LU 分解を用いることが多いが,並列計算の場合には,前処理計算の並列処理手法が必要となる. 本ライブラリルーチンでは並列前処理法として点ヤコビ法,ブロックヤコビ法,加法的シュワルツ法を採用している.
前処理行列Kとして係数行列Aの対角成分のみを用いる.不完全LU分解など他の前処理法と比べると効果は小さいが, 点ヤコビ法は並列に処理され,プロセッサ間通信も不要である.
正方行列を下三角行列と上三角行列の積の形に分解することをLU分解と呼ぶ. LU分解は元の行列が疎行列であっても,分解後の行列の零要素位置に値が入る.これはフィルインと呼ばれる. フィルインが発生することでLU分解後の行列が密行列になるためメモリ使用量が増大する問題がある. PARCELではフィルインを考慮しない不完全LU分解(ILU(0))を前処理として採用することでメモリ使用量の増大を防いでいる.
係数行列をその対角行列,下三角行列,上三角行列を用いて と分解する.この,および対角行列 ()を用いて,前処理行列を の形に決定する.ここでは,対角行列の決定方法として,次の2とおりの条件のどちらを選択する.
条件1 対角成分が一致
条件2 行内の要素和が一致
対角行列の具体的な計算方法は
条件1 対角成分が一致
条件2 行内の要素和が一致
となる.
前処理行列K として係数行列A のブロック対角成分をとり,ブロック内で不完全LU 分解を用いる. 各ブロックの境界を係数行列の各PE への分割と一致させることにより,通信が不要となり,完全に並列に処理できる.
前処理行列K として重なりのあるブロック対角行列を用いる.ブロック内で不完全LU 分解を用いる. 通信は必要となるが,計算は並列に行なわれる. PARCELでは加法的シュワルツ前処理の手法としてBASIC,RESTRICT,INTERPOLATE,NONEを選択できる.
図1 加法的シュワルツ前処理 |
を重なりのある領域まで拡張させたを解き,重なりのある領域のを加算した結果を採用する.
図2 BASIC |
を重なりのある領域まで拡張させたを解き,重なりのある領域のを考慮しない結果を採用する.
図3 RESTRICT |
を重なりのある領域まで拡張させないでを解き,重なりのある領域のを加算した結果を採用する.
図4 INTERPOLATE |
を重なりのある領域まで拡張させないでを解き,重なりのある領域のを考慮しない結果を採用する.
図5 NONE |
一般に疎行列反復解法のプログラムでは,メモリ節約のため,非零要素についての情報(値,行位置,列位置) のみを格納するデータ構造を採用する.ここでは,本ライブラリルーチンで用いている疎行列格納法について説明する.
Compressed Row 形式(圧縮行形式)は,非ゼロ要素を行方向に圧縮した後,列方向にデータを並べて格納する方法である. 非ゼロ要素の列位置・値ともに全体を1 次元配列に格納する.そのため.1次元配列中での各行の開始位置を記録したポインタテーブルを用いる. 並列処理の場合には列方向に行列を分割して格納する.
行列の非ゼロ要素 ::
ポインタテーブル ::
非ゼロ要素の列番号 ::
行列の非ゼロ要素 ::
ポインタテーブル ::
非ゼロ要素の列番号 ::
行列の非ゼロ要素 ::
ポインタテーブル ::
非ゼロ要素の列番号 ::
Diagonal 形式(圧縮対角形式)は,対角行単位に格納する方法である.対角行の値と対角からのオフセットを1次元配列に格納する. 並列処理の場合には列方向に行列を分割して格納する.
対角要素 ::
オフセット ::
対角要素の本数 ::
対角要素 ::
オフセット ::
対角要素の本数 ::
対角要素 ::
オフセット ::
対角要素の本数 ::
PARCELの導入方法について説明する.
展開したアーカイブの内容を図6に示す.
図6 PARCELディレクトリ構成 |
PARCELのコンパイルにはCコンパイラ及びプリプロセッサ対応のFORTRANコンパイラ. MPI,LAPACK,OpenMPライブラリを必要とする. gfortranでコンパイルする場合はコンパイルオプションに-cppを追加することでプリプロセッサが適用される. コンパイルオプションはarch/make_configを編集することで変更する. make_configの記述例を図7に示す. コンパイル方法は展開したアーカイブのトップディレクトリでmakeコマンドを実行する. 再コンパイルする場合はmake cleanを実行後にmakeコマンドを実行する. 各ディレクトリのexampleファイルに実行ファイルが生成されていればコンパイル成功である. 20億次元以上の大規模問題では64ビット整数が必要となる.
図7 make_configの記述例 |
srcディレクトリに生成される*.aファイルをリンクすることで PARCELのサブルーチンを使用することが可能となる.
PARCELに整備されているサブルーチンの直接呼び出しインターフェースについて説明する.
call parcel_dcg( icomm, vecx,vecb, n,gn,nnz,istart, crsA,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret )
共役勾配法(CG法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
call parcel_dcg_dia( icomm, vecx,vecb, n,gn,istart, diaA,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret )
共役勾配法(CG法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
call parcel_ddcg( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,nnz,istart, crsA_hi,crsA_lo,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
共役勾配法(CG法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA_hi(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(上位ビット) |
crsA_lo(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(下位ビット) |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_ddcg_dia( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,istart, diaA_hi,diaA_lo,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
共役勾配法(CG法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA_hi(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(上位ビット) |
diaA_lo(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(下位ビット) |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_dbicgstab( icomm, vecx,vecb, n,gn,nnz,istart, crsA,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret )
安定化双共役勾配法(Bi-CGSTAB法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
call parcel_dbicgstab_dia( icomm, vecx,vecb, n,gn,istart, diaA,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret )
安定化双共役勾配法(Bi-CGSTAB法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
call parcel_ddbicgstab( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,nnz,istart, crsA_hi,crsA_lo,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
安定化双共役勾配法(Bi-CGSTAB法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA_hi(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(上位ビット) |
crsA_lo(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(下位ビット) |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_ddbicgstab_dia( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,istart, diaA_hi,diaA_lo,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
安定化双共役勾配法(Bi-CGSTAB法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA_hi(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(上位ビット) |
diaA_lo(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(下位ビット) |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_dgmres( icomm, vecx,vecb, n,gn,nnz,istart, crsA,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, GMRES_m,GMRES_GSflag )
一般化最小残差法(GMRES(m)法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
GMRES_m | integer | in | リスタートするまでの反復回数 |
GMRES_GSflag | integer | in | 直交化フラグ(1:MGS,2:CGS) |
call parcel_dgmres_dia( icomm, vecx,vecb, n,gn,istart, diaA,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, GMRES_m,GMRES_GSflag )
一般化最小残差法(GMRES(m)法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
GMRES_m | integer | in | リスタートするまでの反復回数 |
GMRES_GSflag | integer | in | 直交化フラグ(1:MGS,2:CGS) |
call parcel_ddgmres( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,nnz,istart, crsA_hi,crsA_lo,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, GMRES_m,GMRES_GSflag, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
一般化最小残差法(GMRES(m)法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA_hi(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(上位ビット) |
crsA_lo(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素(下位ビット) |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
GMRES_m | integer | in | リスタートするまでの反復回数 |
GMRES_GSflag | integer | in | 直交化フラグ(1:MGS,2:CGS) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_ddgmres_dia( icomm, vecx_hi,vecx_lo, vecb_hi,vecb_lo, n,gn,istart, diaA_hi,diaA_lo,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, PRECISION_A, PRECISION_b, PRECISION_x, PRECISION_PRECONDITION )
一般化最小残差法(GMRES(m)法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx_hi(n) | double precision | in/out | in:反復の初期値(上位ビット),out:解ベクトル(上位ビット) |
vecx_lo(n) | double precision | in/out | in:反復の初期値(下位ビット),out:解ベクトル(下位ビット) |
vecb_hi(n) | double precision | in | 右辺ベクトル(上位ビット) |
vecb_lo(n) | double precision | in | 右辺ベクトル(下位ビット) |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA_hi(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(上位ビット) |
diaA_lo(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素(下位ビット) |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
GMRES_m | integer | in | リスタートするまでの反復回数 |
GMRES_GSflag | integer | in | 直交化フラグ(1:MGS,2:CGS) |
PRECISION_A | integer | in | 行列の非零要素の精度(1:倍精度,2:4倍精度) |
PRECISION_b | integer | in | 右辺ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_x | integer | in | 反復の初期値,解ベクトルの精度(1:倍精度,2:4倍精度) |
PRECISION_PRECONDITION | integer | in | 前処理行列の精度(1:倍精度,2:4倍精度) |
call parcel_dcbcg( icomm, vecx,vecb, n,gn,nnz,istart, crsA,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, CBCG_kstep,CBCG_Eigenflag,power_method_itrmax, CAARNOLDI_sstep,CAARNOLDI_tstep, CAARNOLDI_basis,CAARNOLDI_QRflag )
チェビシェフ基底共役勾配法(CBCG法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
CBCG_kstep | integer | in | CBCG法の省通信ステップ |
CBCG_Eigenflag | integer | in | 固有値計算フラグ(1:べき乗法,2:CA-ARNOLDI) |
power_method_itrmax | integer | in | べき乗法の反復回数 |
CAARNOLDI_sstep | integer | in | CA-ARNOLDI法の省通信ステップ数 |
CAARNOLDI_tstep | integer | in | CA-ARNOLDI法のリスタートステップ数 |
CAARNOLDI_basis | integer | in | CA-ARNOLDI法の基底生成フラグ(0:単基底,1:ニュートン基底) |
CAARNOLDI_QRflag | integer | in | CA-ARNOLDI法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
call parcel_dcbcg_dia( icomm, vecx,vecb, n,gn,istart, diaA,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, CBCG_kstep,CBCG_Eigenflag,power_method_itrmax, CAARNOLDI_sstep,CAARNOLDI_tstep, CAARNOLDI_basis,CAARNOLDI_QRflag )
チェビシェフ基底共役勾配法(CBCG法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
CBCG_kstep | integer | in | CBCG法の省通信ステップ |
CBCG_Eigenflag | integer | in | 固有値計算フラグ(1:べき乗法,2:CA-ARNOLDI) |
power_method_itrmax | integer | in | べき乗法の最大反復回数 |
CAARNOLDI_sstep | integer | in | CA-ARNOLDI法の省通信ステップ |
CAARNOLDI_tstep | integer | in | CA-ARNOLDI法のリスタートステップ |
CAARNOLDI_basis | integer | in | CA-ARNOLDI法の基底生成フラグ(0:単基底,1:ニュートン基底) |
CAARNOLDI_QRflag | integer | in | CA-ARNOLDI法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
call parcel_dcagmres( icomm, vecx,vecb, n,gn,nnz,istart, crsA,crsRow_ptr,crsCol, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, CAGMRES_sstep, CAGMRES_tstep, CAGMRES_basis, CAGMRES_QRflag )
通信削減型一般化最小残差法(CA-GMRES(s,t)法)により,連立一次方程式Ax=bを解く.非零要素の値をCompressed Row形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
CAGMRES_sstep | integer | in | CA-GMRES法の省通信ステップ数 |
CAGMRES_tstep | integer | in | CA-GMRES法のリスタートステップ数 |
CAGMRES_basis | integer | in | CA-GMRES法の基底生成フラグ(0:単基底,1:ニュートン基底) |
CAGMRES_QRflag | integer | in | CA-GMRES法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
call parcel_dcagmres_dia( icomm, vecx,vecb, n,gn,istart, diaA,offset,nnd, ipreflag,ILU_method,addL,iflagAS, itrmax,rtolmax, reshistory, iovlflag,iret, CAGMRES_sstep, CAGMRES_tstep, CAGMRES_basis, CAGMRES_QRflag )
通信削減型一般化最小残差法(CA-GMRES(s,t)法)により,連立一次方程式Ax=bを解く.非零要素の値をDiagonal形式で格納.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in/out | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
iflagAS | integer | in | 加法シュワルツ指定フラグ(1:BASIC,2:RESTRICT,3:INTERPOLATE,4:NONE) |
itrmax | integer | in/out | in:最大反復回数,out:反復回数 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
reshistory(itrmax) | double precision | out | 残差履歴 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
iret | integer | out | エラーコード(0:正常終了) |
CAGMRES_sstep | integer | in | CA-GMRES法の省通信ステップ |
CAGMRES_tstep | integer | in | CA-GMRES法のリスタートステップ |
CAGMRES_basis | integer | in | CA-GMRES法の基底生成フラグ(0:単基底,1:ニュートン基底) |
CAGMRES_QRflag | integer | in | CA-GMRES法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
PARCELに整備されているサブルーチンの間接呼び出しインターフェースについて説明する.
call PARCEL_KSP_Default_Setting(method)
PARCEL構造体に初期設定を適用.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
call set_KSP_CRS(icomm,vecx,vecb,n,gn,nnz,istart,crsA,crsRow_ptr,crsCol,itrmax,ipreflag,ILU_method,addL,method)
Compressed Row形式の行列を設定.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
nnz | integer*4 / integer*8 | in | 各プロセスが担当する係数行列の非零要素数 |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
crsA(nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
crsRow_ptr(n+1) | integer*4 / integer*8 | in | Compressed Row 形式のポインタテーブル |
crsCol(nnz) | integer*4 / integer*8 | in | Compressed Row 形式に格納された分割行列の非零要素の列番号 |
itrmax | integer | in | in:最大反復回数,out:反復回数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
method | type(KSP) | in/out | PARCEL 構造体 |
call set_KSP_DIA(icomm,vecx,vecb,n,gn,istart,diaA,offset,nnd,itrmax,ipreflag,ILU_method,addL,method)
Diagonal形式の行列を設定.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
icomm | integer | in | MPIコミュニケータ |
vecx(n) | double precision | in | in:反復の初期値,out:解ベクトル |
vecb(n) | double precision | in | 右辺ベクトル |
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
gn | integer*4 / integer*8 | in | 全体のベクトルの大きさ |
istart | integer*4 / integer*8 | in | 各プロセスが担当する行列の開始行 |
diaA(n*nnd) | double precision | in | Diagonal形式に格納された行列の非零要素 |
offset(nnd) | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角からのオフセット |
nnd | integer*4 / integer*8 | in | Diagonal形式に格納された行列の対角要素の本数 |
itrmax | integer | in | in:最大反復回数,out:反復回数 |
ipreflag | integer | in | 前処理指定フラグ(0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(0:ILU(0),1:D-ILU(対角成分が一致),2:D-ILU(行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅 |
method | type(KSP) | in/out | PARCEL 構造体 |
call set_vecx_KSP(n,x,method)
反復の初期値を設定.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
x(n) | double precision | in | 反復の初期値 |
method | type(KSP) | in/out | PARCEL 構造体 |
call set_vecb_KSP(n,b,method)
右辺ベクトルbを設定.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
b(n) | double precision | in | 右辺ベクトル |
method | type(KSP) | in/out | PARCEL 構造体 |
call get_vecx_KSP(n,x,method)
解ベクトルを取得.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
n | integer*4 / integer*8 | in | 各プロセスが担当するベクトルの大きさ |
x(n) | double precision | out | 解ベクトル |
method | type(KSP) | in | PARCEL 構造体 |
call get_reshistory_KSP(niter,reshistory, method)
収束履歴を取得する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
niter | integer | out | 収束までの反復回数 |
reshistory(method%itrmax) | double precision | out | 収束履歴 |
method | type(KSP) | in | PARCEL 構造体 |
call reset_CRS_Mat(itrmax,ipreflag,ILU_method,addL,crsA,method)
行列の非零要素位置を変えずにCompressed Row形式の行列値,前処理方法,最大反復回数を変更する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
itrmax | integer | in | 最大反復回数(-1:設定済みの値を引き継ぐ) |
ipreflag | integer | in | 前処理指定フラグ(-1:設定済みの値を引き継ぐ,0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(-1:設定済みの値を引き継ぐ,1:対角成分が一致,2:行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅(-1:設定済みの値を引き継ぐ) |
crsA(method%nnz) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
method | type(KSP) | in/out | PARCEL 構造体 |
call reset_DIA_Mat(itrmax,ipreflag,ILU_method,addL,diaMat,method)
行列の非零要素位置を変えずにDiagonal形式の行列値,前処理方法,最大反復回数を変更する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
itrmax | integer | in | 最大反復回数(-1:設定済みの値を引き継ぐ) |
ipreflag | integer | in | 前処理指定フラグ(-1:設定済みの値を引き継ぐ,0:無し,1:ポイントヤコビ,2:ブロックヤコビ,3:加法シュワルツ) |
ILU_method | integer | in | 不完全LU分解指定フラグ(-1:設定済みの値を引き継ぐ,1:対角成分が一致,2:行内の要素和が一致) |
addL | integer | in | 加法シュワルツの重なり幅(-1:設定済みの値を引き継ぐ) |
diaMat(method%n*method%nnd) | double precision | in | Compressed Row 形式に格納された行列の非零要素 |
method | type(KSP) | in/out | PARCEL 構造体 |
call set_KSP_solver(method,solver)
クリロフ部分空間法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
solver | integer | in | クリロフ部分空間法フラグ(1:CG法,2:Bi-CGSTAB法,3:GMRES法,4:CAGMRES法,5:CBCG法) |
call set_KSP_abstolmax(method,abstolmax)
クリロフ部分空間法の絶対残差ノルムによる収束判定条件を設定する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
abstolmax | double precision | in | 収束判定条件(絶対残差ノルム) |
call set_KSP_rtolmax(method,rtolmax)
クリロフ部分空間法の相対残差ノルムによる収束判定条件を設定する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
rtolmax | double precision | in | 収束判定条件(相対残差ノルム) |
call set_KSP_iovlflag(method,iovlflag)
疎行列ベクトル積の通隠蔽手法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
iovlflag | integer | in | 通信隠蔽フラグ(0:無し,1:全プロセス,2:各プロセス) |
call set_KSP_GMRES_m(method,GMRES_m)
GMRES(m)のリスタート数mを選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
GMRES_m | integer | in | リスタートするまでの反復回数 |
call set_KSP_GMRES_GSflag(method,GMRES_GSflag)
GMRES(m)の直行化手法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
GMRES_GSflag | integer | in | 直交化フラグ(1:MGS,2:CGS) |
call set_KSP_CAGMRES_sstep(method,CAGMRES_sstep)
CA-GMRES法の省通信ステップ数(s)を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAGMRES_sstep | integer | in | CA-GMRES法の省通信ステップ数 |
call set_KSP_CAGMRES_tstep(method,CAGMRES_tstep)
CA-GMRES法のリスタートステップ数(t)を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAGMRES_tstep | integer | in | CA-GMRES法のリスタートステップ数 |
call set_KSP_CAGMRES_basis(method,CAGMRES_basis)
CA-GMRES法の基底生成方法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAGMRES_basis | integer | in | CA-GMRES法の基底生成フラグ(0:単基底,1:ニュートン基底) |
call set_KSP_CAGMRES_QRflag(method,CAGMRES_QRflag)
CA-GMRES法のQR分解方法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAGMRES_QRflag | integer | in | CA-GMRES法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
call set_KSP_CBCG_kstep(method,CBCG_kstep)
CBCG法の省通信ステップ数kを選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CBCG_kstep | integer | in | CBCG法の省通信ステップ |
call set_KSP_CBCG_Eigenflag(method,CBCG_Eigenflag)
CBCG法の固有値計算方法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CBCG_Eigenflag | integer | in | 固有値計算フラグ(1:べき乗法,2:CA-ARNOLDI) |
call set_KSP_CAARNOLDI_sstep(method,CAARNOLDI_sstep)
CA-ARNOLDI法の省通信ステップ数(s)を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAARNOLDI_sstep | integer | in | CA-ARNOLDI法の省通信ステップ数 |
call set_KSP_CAARNOLDI_tstep(method,CAARNOLDI_tstep)
CA-ARNOLDI法のリスタートステップ(t)数を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAARNOLDI_tstep | integer | in | CA-ARNOLDI法のリスタートステップ |
call set_KSP_CAARNOLDI_basis(method,CAARNOLDI_basis)
CA-ARNOLDI法の基底生成方法を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAARNOLDI_basis | integer | in | CA-ARNOLDI法の基底生成フラグ(0:単基底,1:ニュートン基底) |
call set_KSP_CAARNOLDI_QRflag(method,CAARNOLDI_QRflag)
CA-ARNOLDI法のQR分解を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
CAARNOLDI_QRflag | integer | in | CA-ARNOLDI法のQR分解フラグ(1:MGS,2:CGS,3:TSQR,4:CholeskyQR,5:CholeskyQR2) |
call set_KSP_power_method_itrmax(method,power_method_itrmax)
べき乗法の反復回数を選択する.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | out | PARCEL 構造体 |
power_method_itrmax | integer | in | べき乗法の反復回数 |
call START_KSP_SOLVER(method)
連立一次方程式Ax=bの求解.
引数名(次元) | 型 | 入力/出力 | 説明 |
---|---|---|---|
method | type(KSP) | in/out | PARCEL 構造体 |
PARCELのサブルーチン使用方法についてCG法のFORTRANサンプルコードを例にして説明する. Cプログラムにおける利用法の詳細についてはアーカイブのSPARCE/example_Cに格納されているCサンプルコードを参照のこと. Cの場合はFORTRAN用のMPIコミュニケータを用意すればFORTRANと同様の方法でPARCELを利用可能である. CサンプルコードではMPIライブラリの関数MPI_Comm_c2fを用いてCのMPIコミュニケータをFORTRAN用に変換している. 間接呼び出しインターフェースのCサンプルコードでは,FORTRAN側でPARCEL構造体を確保してからC側でそれを制御して利用している.
Compressed Row形式のサンプルコードを図8および図9に示す. make_matrix_CRSはPARCELの行列分割に対応したCompressed Row形式の行列を生成する任意のサブルーチンとする. 図8 は直接呼び出しインターフェースを用いた例,図9は間接呼び出しインターフェースを利用した例を示す. 図9ではkrylov,krylov_kernelのモジュールを利用するためSPARSE/srcをコンパイル時にインクルードする必要がある. 非零要素の位置を変更せずに値のみを変更する場合, 図8のコードでは再度parcel_dcgを呼ぶ必要があるため通信相手のリスト作成等の初期化コストが必要となる. 図9のコードではreset_CRS_Matを利用することで非零要素位置の変更がない場合は初期化コストは不要となる. 同様に右辺ベクトル及び反復の初期値を変更する場合も図8ではparcel_dcgを呼ぶため初期化のコストが必要となるが, 図9のコードではset_vecb_KSP,set_vecx_KSPを利用することで初期化コストは不要となる.
program main use mpi implicit none integer n,gn,nnz,istart real*8,allocatable :: crsA(:) integer,allocatable :: crsRow_ptr(:),crsCol(:) real*8,allocatable :: vecx(:) real*8,allocatable :: vecb(:) integer itrmax real*8 rtolmax real*8, allocatable :: reshistory(:) integer ipreflag integer ILU_method integer iflagAS integer addL integer iovlflag integer iret integer ierr call MPI_Init(ierr) call make_matrix_CRS(n,gn,nnz,istart,crsA,crsRow_ptr,crsCol) allocate(vecx(n)) allocate(vecb(n)) allocate(reshistory(itrmax)) ipreflag=0 ILU_method=1 addL=0 iflagAS=1 itrmax=100 rtolmax=1.0d-8 iovlflag=0 vecb=1.0d0 vecx=1.0d0 call parcel_dcg( & MPI_COMM_WORLD, & vecx,vecb, & n,gn,nnz,istart, & crsA,crsRow_ptr,crsCol, & ipreflag,ILU_method,addL,iflagAS, & itrmax,rtolmax, & reshistory, & iovlflag,iret & ) call MPI_Finalize(ierr) deallocate(vecx) deallocate(vecb) deallocate(reshistory) end program main
図8 Compressed Row形式 サンプルコード(直接呼び出しインターフェース) |
program main use mpi use krylov use krylov_kernel implicit none integer n,gn,nnz,istart real*8,allocatable :: crsA(:) integer,allocatable :: crsRow_ptr(:),crsCol(:) real*8,allocatable :: vecx(:) real*8,allocatable :: vecb(:) integer itrmax real*8 rtolmax real*8, allocatable :: reshistory(:) integer ipreflag integer ILU_method integer iflagAS integer addL integer iovlflag integer iret integer ierr type(KSP) :: method call MPI_Init(ierr) call make_matrix_CRS(n,gn,nnz,istart,crsA,crsRow_ptr,crsCol) allocate(vecx(n)) allocate(vecb(n)) allocate(reshistory(itrmax)) ipreflag=0 ILU_method=1 addL=0 iflagAS=1 itrmax=100 rtolmax=1.0d-8 iovlflag=0 vecb=1.0d0 vecx=1.0d0 call PARCEL_KSP_Default_Setting(method) call set_KSP_CRS(& MPI_COMM_WORLD, & vecx,vecb, & n,gn,nnz,istart, & crsA,crsRow_ptr,crsCol, & itrmax,ipreflag,ILU_method,addL,& method) call set_KSP_solver(method,1) call START_KSP_SOLVER(method) call get_vecx_KSP(n,vecx,method) call free_KSP(method) call MPI_Finalize(ierr) deallocate(vecx) deallocate(vecb) deallocate(reshistory) end program main
図9 Compressed Row形式 サンプルコード(間接呼び出しインターフェース) |
Diagonal形式のサンプルコードを図10に示す. make_matrix_DIAはPARCELのブロック分割に対応したDiagonal形式の行列を生成する任意のサブルーチンとする. 図10は直接呼び出しインターフェースを用いた例, 図11は間接呼び出しインターフェースを利用した例を示す. 図11の形式ではkrylov,krylov_kernelのモジュールを利用するためSPARSE/srcをコンパイル時にインクルードする必要がある. 非零要素の位置を変更せずに値のみを変更する場合, 図10のコードでは再度parcel_dcg_diaを呼ぶ必要があるため通信相手のリスト作成等の初期化コストが必要となる. 図11のコードではreset_DIA_Matを利用することで非零要素位置の変更がない場合は初期化コストは不要となる. 同様に右辺ベクトル及び反復の初期値を変更する場合も図10ではparcel_dcg_diaを呼ぶため初期化のコストが必要となるが, 図11のコードではset_vecb_KSP,set_vecx_KSPを利用することで初期化コストは不要となる.
program main use mpi implicit none integer n,gn,nnd,istart real*8,allocatable :: diaA(:) integer,allocatable :: offset(:) real*8,allocatable :: vecx(:) real*8,allocatable :: vecb(:) integer itrmax real*8 rtolmax real*8, allocatable :: reshistory(:) integer ipreflag integer ILU_method integer iflagAS integer addL integer iovlflag integer iret integer ierr call MPI_Init(ierr) call make_matrix_DIA(n,gn,nnd,istart,diaA,offset) allocate(vecx(n)) allocate(vecb(n)) allocate(reshistory(itrmax)) ipreflag=0 ILU_method=1 addL=0 iflagAS=1 itrmax=100 rtolmax=1.0d-8 iovlflag=0 vecb=1.0d0 vecx=1.0d0 call parcel_dcg_dia(& MPI_COMM_WORLD, & vecx,vecb,& n,gn,istart,& diaA,offset,nnd,& ipreflag,ILU_method,addL,iflagAS,& itrmax,rtolmax,& reshistory,& iovlflag,iret& ) call MPI_Finalize(ierr) deallocate(reshistory) end program main
図10 Diagonal形式 サンプルコード(直接呼び出しインターフェース) |
program main use mpi use krylov use krylov_kernel implicit none integer n,gn,nnd,istart real*8,allocatable :: diaA(:) integer,allocatable :: offset(:) real*8,allocatable :: vecx(:) real*8,allocatable :: vecb(:) integer itrmax real*8 rtolmax real*8, allocatable :: reshistory(:) integer ipreflag integer ILU_method integer iflagAS integer addL integer iovlflag integer iret integer ierr call MPI_Init(ierr) call make_matrix_DIA(n,gn,nnd,istart,diaA,offset) allocate(vecx(n)) allocate(vecb(n)) allocate(reshistory(itrmax)) ipreflag=0 ILU_method=1 addL=0 iflagAS=1 itrmax=100 rtolmax=1.0d-8 iovlflag=0 vecb=1.0d0 vecx=1.0d0 call PARCEL_KSP_Default_Setting(method) call set_KSP_DIA(& MPI_COMM_WORLD, & vecx,vecb, & n,gn,istart, & diaA,offset,nnd,& itrmax,ipreflag,ILU_method,addL,& method) call set_KSP_solver(method,1) call START_KSP_SOLVER(method) call get_vecx_KSP(n,vecx,method) call free_KSP(method) call MPI_Finalize(ierr) deallocate(vecx) deallocate(vecb) deallocate(reshistory) end program main
図11 Diagonal形式 サンプルコード(間接呼び出しインターフェース) |