/////////////////////////////////////////////////////////////////////////////// // positive definite band // Symmetric / Hermitian Positive Definite Band Matrix private bool DoublePositiveDefiniteBandSolve(out double[] X, DoubleSparseMatrix A, double[] B) { X = null; bool isSymmetric = A.IsSymmetric(); System.Diagnostics.Debug.Assert(isSymmetric); if (!isSymmetric) { return(false); } bool success = false; IvyFEM.Lapack.DoubleSymmetricBandMatrix pbA = (IvyFEM.Lapack.DoubleSymmetricBandMatrix)A; A = null; int xRow; int xCol; int ret = IvyFEM.Lapack.Functions.dpbsv(out X, out xRow, out xCol, pbA.Buffer, pbA.RowLength, pbA.ColumnLength, pbA.SuperdiaLength, B, B.Length, 1); System.Diagnostics.Debug.Assert(ret == 0); success = (ret == 0); return(success); }
public static bool DoubleSolveBiCGSTABWithPivoting(out double[] X, DoubleSparseMatrix A, double[] B, int fillinLevel, double convRatioTolerance) { int t; t = System.Environment.TickCount; DoubleSparseMatrix LU; int[] pivot; DoubleCalcILUWithPivoting(out LU, out pivot, A, fillinLevel); DoubleSparseMatrix pivotingA = new DoubleSparseMatrix(A.RowLength, A.ColumnLength); for (int i = 0; i < pivotingA.RowLength; i++) { pivotingA.RowColIndexValues[i] = A.RowColIndexValues[pivot[i]]; } double[] pivotingB = new double[B.Length]; for (int i = 0; i < B.Length; i++) { pivotingB[i] = B[pivot[i]]; } A = null; B = null; System.Diagnostics.Debug.WriteLine("DoubleSolveBiCGStab 1: t= " + (System.Environment.TickCount - t)); t = System.Environment.TickCount; bool success = IvyFEM.Linear.Functions.DoubleSolvePreconditionedBiCGSTAB(out X, pivotingA, pivotingB, LU, convRatioTolerance); System.Diagnostics.Debug.WriteLine("DoubleSolveBiCGSTAB 2: t= " + (System.Environment.TickCount - t)); return(success); }
private bool DoubleSolveBiCGSTABWithPivoting(out double[] X, DoubleSparseMatrix A, double[] B) { int t; X = null; //t = System.Environment.TickCount; //bool success = IvyFEM.Linear.Functions.DoubleSolveBiCGSTABWithPivoting(out X, A, B, ILUFillinLevel, // ConvRatioTolerance); //System.Diagnostics.Debug.Assert(success); //System.Diagnostics.Debug.WriteLine(" DoubleSolveCGWithPivoting t = " + (System.Environment.TickCount - t)); // Nativeを使う bool success = false; { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; int[] APtrs; int[] AIndexs; double[] AValues; t = System.Environment.TickCount; A.GetCSR(out APtrs, out AIndexs, out AValues); System.Diagnostics.Debug.WriteLine(" GetCSR t = " + (System.Environment.TickCount - t)); t = System.Environment.TickCount; success = IvyFEM.Native.Functions.DoubleSolveBiCGSTABWithPivoting( out X, n, APtrs, AIndexs, AValues, B, ILUFillinLevel, ConvRatioTolerance); System.Diagnostics.Debug.Assert(success); System.Diagnostics.Debug.WriteLine(" DoubleSolveCGWithPivoting t = " + (System.Environment.TickCount - t)); } return(success); }
public bool DoubleSolve(out double[] X, DoubleSparseMatrix A, double[] B) { bool success = false; X = null; switch (Method) { case LapackEquationSolverMethod.Dense: success = DoubleDenseSolve(out X, A, B); break; case LapackEquationSolverMethod.Default: case LapackEquationSolverMethod.Band: success = DoubleBandSolve(out X, A, B); break; case LapackEquationSolverMethod.PositiveDefiniteBand: success = DoublePositiveDefiniteBandSolve(out X, A, B); break; default: throw new NotImplementedException(); //break; } return(success); }
///////////////////////////////////////////////////////////////////////////////// // double // LU : Lii = 1 Lij (i=1...n -1 j=0...i-1) Uij (i=0...n-1 j=i...n-1) // Note: M = L * U public static DoubleSparseMatrix DoubleCalcILU(DoubleSparseMatrix A, int fillinLevel) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; DoubleSparseMatrix LU = new DoubleSparseMatrix(A); int[,] level = new int[n, n]; for (int row = 0; row < n; row++) { for (int col = 0; col < n; col++) { level[row, col] = (Math.Abs(A[row, col]) >= IvyFEM.Constants.PrecisionLowerLimit) ? 0 : (fillinLevel + 1); } } for (int i = 1; i < n; i++) { for (int k = 0; k <= (i - 1); k++) { //if (!LU.RowColIndexValues[i].ContainsKey(k)) //{ // continue; //} if (level[i, k] > fillinLevel) { continue; } LU[i, k] /= LU[k, k]; double LUik = LU[i, k]; foreach (var pair in LU.RowColIndexValues[k]) { int j = pair.Key; double LUkj = pair.Value; if (j >= k + 1 && j < n) { // } else { continue; } level[i, j] = Math.Min(level[i, j], level[i, k] + level[k, j] + 1); if (level[i, j] <= fillinLevel) { LU[i, j] -= LUik * LUkj; } } } } return(LU); }
public void Copy(DoubleSparseMatrix src) { Resize(src.RowLength, src.ColumnLength); for (int row = 0; row < src.RowLength; row++) { var srcColIndexValues = src.RowColIndexValues[row]; foreach (var pair in srcColIndexValues) { int col = pair.Key; double value = pair.Value; RowColIndexValues[row].Add(col, value); } } }
/////////////////////////////////////////////////////////////////////////////// // dense private bool DoubleDenseSolve(out double[] X, DoubleSparseMatrix A, double[] B) { bool success = false; IvyFEM.Lapack.DoubleMatrix denseA = (IvyFEM.Lapack.DoubleMatrix)A; int xRow; int xCol; int ret = IvyFEM.Lapack.Functions.dgesv(out X, out xRow, out xCol, denseA.Buffer, denseA.RowLength, denseA.ColumnLength, B, B.Length, 1); System.Diagnostics.Debug.Assert(ret == 0); success = (ret == 0); return(success); }
private static bool[,] GetDoubleMatrixNonzeroPattern(DoubleSparseMatrix A) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; bool[,] nonzeroPattern = new bool[n, n]; for (int row = 0; row < n; row++) { foreach (var pair in A.RowColIndexValues[row]) { int col = pair.Key; nonzeroPattern[row, col] = true; } } return(nonzeroPattern); }
public static bool DoubleSolveICCG(out double[] X, DoubleSparseMatrix A, double[] B, double convRatioTolerance) { int t; t = System.Environment.TickCount; DoubleSparseMatrix LU = IvyFEM.Linear.Functions.DoubleCalcIC(A); System.Diagnostics.Debug.WriteLine("DoubleSolveICCG 1: t= " + (System.Environment.TickCount - t)); t = System.Environment.TickCount; bool success = IvyFEM.Linear.Functions.DoubleSolvePreconditionedCG(out X, A, B, LU, convRatioTolerance); System.Diagnostics.Debug.WriteLine("DoubleSolveICCG 2: t= " + (System.Environment.TickCount - t)); return(success); }
public static bool DoubleSolveBiCGSTAB(out double[] X, DoubleSparseMatrix A, double[] B, int fillinLevel, double convRatioTolerance) { int t; t = System.Environment.TickCount; DoubleSparseMatrix LU = IvyFEM.Linear.Functions.DoubleCalcILU(A, fillinLevel); System.Diagnostics.Debug.WriteLine("DoubleSolveBiCGStab 1: t= " + (System.Environment.TickCount - t)); t = System.Environment.TickCount; bool success = IvyFEM.Linear.Functions.DoubleSolvePreconditionedBiCGSTAB(out X, A, B, LU, convRatioTolerance); System.Diagnostics.Debug.WriteLine("DoubleSolveBiCGSTAB 2: t= " + (System.Environment.TickCount - t)); return(success); }
private bool __DoubleBandSolve(out double[] X, DoubleSparseMatrix A, double[] B) { bool success = false; IvyFEM.Lapack.DoubleBandMatrix bandA = (IvyFEM.Lapack.DoubleBandMatrix)A; A = null; int xRow; int xCol; int ret = IvyFEM.Lapack.Functions.dgbsv(out X, out xRow, out xCol, bandA.Buffer, bandA.RowLength, bandA.ColumnLength, bandA.SubdiaLength, bandA.SuperdiaLength, B, B.Length, 1); System.Diagnostics.Debug.Assert(ret == 0); success = (ret == 0); return(success); }
public bool DoubleSolve(out double[] X, DoubleSparseMatrix A, double[] B) { bool success = false; X = null; switch (Method) { case IvyFEMEquationSolverMethod.Default: case IvyFEMEquationSolverMethod.NoPreconCG: success = DoubleSolveNoPreconCG(out X, A, B); break; case IvyFEMEquationSolverMethod.CG: success = DoubleSolveCG(out X, A, B); break; case IvyFEMEquationSolverMethod.CGWithPivoting: success = DoubleSolveCGWithPivoting(out X, A, B); break; case IvyFEMEquationSolverMethod.ICCG: success = DoubleSolveICCG(out X, A, B); break; case IvyFEMEquationSolverMethod.NoPreconBiCGSTAB: success = DoubleSolveNoPreconBiCGSTAB(out X, A, B); break; case IvyFEMEquationSolverMethod.BiCGSTAB: success = DoubleSolveBiCGSTAB(out X, A, B); break; case IvyFEMEquationSolverMethod.BiCGSTABWithPivoting: success = DoubleSolveBiCGSTABWithPivoting(out X, A, B); break; default: throw new NotImplementedException(); //break; } return(success); }
private bool DoubleSolveICCG(out double[] X, DoubleSparseMatrix A, double[] B) { int t; X = null; //t = System.Environment.TickCount; //bool isSymmetric = A.IsSymmetric(); //System.Diagnostics.Debug.Assert(isSymmetric); //System.Diagnostics.Debug.WriteLine(" IsSymmetric t = " + (System.Environment.TickCount - t)); //if (!isSymmetric) //{ // return false; //} //t = System.Environment.TickCount; //bool success = IvyFEM.Linear.Functions.DoubleSolveICCG(out X, A, B, // ConvRatioTolerance); //System.Diagnostics.Debug.Assert(success); //System.Diagnostics.Debug.WriteLine(" DoubleSolveICCG t = " + (System.Environment.TickCount - t)); // Nativeを使う bool success = false; { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; int[] APtrs; int[] AIndexs; double[] AValues; t = System.Environment.TickCount; A.GetCSR(out APtrs, out AIndexs, out AValues); System.Diagnostics.Debug.WriteLine(" GetCSR t = " + (System.Environment.TickCount - t)); t = System.Environment.TickCount; success = IvyFEM.Native.Functions.DoubleSolveICCG( out X, n, APtrs, AIndexs, AValues, B, ConvRatioTolerance); System.Diagnostics.Debug.Assert(success); System.Diagnostics.Debug.WriteLine(" DoubleSolveICCG t = " + (System.Environment.TickCount - t)); } return(success); }
/////////////////////////////////////////////////////////////////////////////// // band private bool DoubleBandSolve(out double[] X, DoubleSparseMatrix A, double[] B) { X = null; DoubleSparseMatrix AToSolve; double[] BToSolve; int[] indexs; if (IsOrderingToBandMatrix) { bool successOrder = IvyFEM.Linear.Utils.OrderToDoubleBandMatrix( out AToSolve, out BToSolve, out indexs, A, B); System.Diagnostics.Debug.Assert(successOrder); if (!successOrder) { return(false); } } else { AToSolve = A; BToSolve = B; indexs = null; } double[] XToSolve; bool success = __DoubleBandSolve(out XToSolve, AToSolve, BToSolve); if (IsOrderingToBandMatrix) { X = new double[XToSolve.Length]; for (int row = 0; row < XToSolve.Length; row++) { X[indexs[row]] = XToSolve[row]; } } else { X = XToSolve; } return(success); }
public static void DoubleSolveLU(out double[] X, DoubleSparseMatrix LU, double[] B) { System.Diagnostics.Debug.Assert(LU.RowLength == LU.ColumnLength); int n = LU.RowLength; X = new double[n]; B.CopyTo(X, 0); // Ly = b for (int row = 1; row < n; row++) { foreach (var pair in LU.RowColIndexValues[row]) { int col = pair.Key; double value = pair.Value; if (col >= 0 && col < row) { X[row] -= value * X[col]; // LU[row, col] } } } // Ux = y for (int row = n - 2; row >= 0; row--) { foreach (var pair in LU.RowColIndexValues[row]) { int col = pair.Key; double value = pair.Value; if (col >= row + 1 && col < n) { X[row] -= value * X[col]; // LU[row, col] } } X[row] /= LU[row, row]; } }
/* * public void Transpose() * { * throw new NotImplementedException(); * } */ public static DoubleSparseMatrix operator *(DoubleSparseMatrix A, DoubleSparseMatrix B) { int aRow = A.RowLength; int aCol = A.ColumnLength; int bRow = B.RowLength; int bCol = B.ColumnLength; if (aCol != bRow) { throw new ArgumentException("Mismatched size: aCol != bRow(" + aCol + " != " + bRow + ")"); } int cRow = aRow; int cCol = bCol; DoubleSparseMatrix C = new DoubleSparseMatrix(cRow, cCol); for (int i = 0; i < aRow; i++) { var colIndexValues1 = A.RowColIndexValues[i]; for (int j = 0; j < bCol; j++) { foreach (var pair1 in colIndexValues1) { int k = pair1.Key; double value1 = pair1.Value; // Aik if (B.RowColIndexValues[k].ContainsKey(j)) { double value2 = B.RowColIndexValues[k][j]; // Bkj C[i, j] += value1 * value2; } } } } return(C); }
public bool DoubleSolve(out double[] X, DoubleSparseMatrix A, double[] B) { bool success = false; using (IvyFEM.Lis.LisInitializer LisInitializer = new IvyFEM.Lis.LisInitializer()) using (IvyFEM.Lis.LisMatrix lisA = (IvyFEM.Lis.LisMatrix)A) using (IvyFEM.Lis.LisVector lisB = (IvyFEM.Lis.LisVector)B) using (IvyFEM.Lis.LisVector lisX = new IvyFEM.Lis.LisVector()) using (IvyFEM.Lis.LisSolver lisSolver = new IvyFEM.Lis.LisSolver()) { int ret; int n = B.Length; lisX.SetSize(0, n); A = null; ret = lisSolver.SetOption(GetOption()); System.Diagnostics.Debug.Assert(ret == 0); ret = lisSolver.SetOptionC(); System.Diagnostics.Debug.Assert(ret == 0); ret = lisSolver.Solve(lisA, lisB, lisX); System.Diagnostics.Debug.Assert(ret == 0); success = (ret == 0); int iter; ret = lisSolver.GetIter(out iter); System.Diagnostics.Debug.Assert(ret == 0); System.Diagnostics.Debug.WriteLine("Lis Solve iter = " + iter); System.Numerics.Complex[] complexX = new System.Numerics.Complex[n]; ret = lisX.GetValues(0, n, complexX); System.Diagnostics.Debug.Assert(ret == 0); X = new double[n]; for (int i = 0; i < n; i++) { X[i] = complexX[i].Real; } } return(success); }
public static bool DoubleSolveNoPreconCG( out double[] X, DoubleSparseMatrix A, double[] B, double convRatioTolerance) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; System.Diagnostics.Debug.Assert(B.Length == n); double convRatio = convRatioTolerance; double tolerance = convRatio; const int maxIter = IvyFEM.Linear.Constants.MaxIter; double[] r = new double[n]; double[] x = new double[n]; double[] p = new double[n]; int iter = 0; B.CopyTo(r, 0); double rr; double sqInvNorm0; { double sqNorm0 = IvyFEM.Lapack.Functions.ddot(r, r); if (sqNorm0 < IvyFEM.Constants.PrecisionLowerLimit) { convRatio = 0; X = x; System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); return(true); } rr = sqNorm0; sqInvNorm0 = 1.0 / sqNorm0; } r.CopyTo(p, 0); for (iter = 0; iter < maxIter; iter++) { double[] Ap = A * p; double alpha; { double pAp = IvyFEM.Lapack.Functions.ddot(p, Ap); alpha = rr / pAp; } r = IvyFEM.Lapack.Functions.daxpy(-alpha, Ap, r); x = IvyFEM.Lapack.Functions.daxpy(alpha, p, x); double rrPrev = rr; { double sqNorm = IvyFEM.Lapack.Functions.ddot(r, r); if (sqNorm * sqInvNorm0 < tolerance * tolerance) { convRatio = Math.Sqrt(sqNorm * sqInvNorm0); X = x; System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); return(true); } rr = sqNorm; } double beta = rr / rrPrev; p = IvyFEM.Lapack.Functions.daxpy(beta, p, r); } { double sqNormRes = IvyFEM.Lapack.Functions.ddot(r, r); convRatio = Math.Sqrt(sqNormRes * sqInvNorm0); System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); X = x; } System.Diagnostics.Debug.WriteLine("Not converged"); return(false); }
public static void DoubleCalcILUWithPivoting(out DoubleSparseMatrix LU, out int[] pivot, DoubleSparseMatrix A, int fillinLevel) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; LU = new DoubleSparseMatrix(A); pivot = new int[n]; for (int row = 0; row < n; row++) { pivot[row] = row; } int[,] level = new int[n, n]; for (int row = 0; row < n; row++) { for (int col = 0; col < n; col++) { level[row, col] = (Math.Abs(A[row, col]) >= IvyFEM.Constants.PrecisionLowerLimit) ? 0 : (fillinLevel + 1); } } for (int k = 0; k < (n - 1); k++) { int p = k; { double max = Math.Abs(LU[k, k]); for (int i = k + 1; i < n; i++) { double abs = Math.Abs(LU[i, k]); if (abs > max) { max = abs; p = i; } } } if (k != p) { { int tmp = pivot[k]; pivot[k] = pivot[p]; pivot[p] = tmp; } { var tmp = LU.RowColIndexValues[k]; LU.RowColIndexValues[k] = LU.RowColIndexValues[p]; LU.RowColIndexValues[p] = tmp; } for (int j = 0; j < n; j++) { int tmp = level[k, j]; level[k, j] = level[p, j]; level[p, j] = tmp; } } for (int i = k + 1; i < n; i++) { if (level[i, k] > fillinLevel) { continue; } System.Diagnostics.Debug.Assert(LU.RowColIndexValues[k].ContainsKey(k)); LU[i, k] /= LU[k, k]; double LUik = LU[i, k]; foreach (var pair in LU.RowColIndexValues[k]) { int j = pair.Key; double LUkj = pair.Value; if (j >= k + 1 && j < n) { } else { continue; } level[i, j] = Math.Min(level[i, j], level[i, k] + level[k, j] + 1); if (level[i, j] <= fillinLevel) { LU[i, j] -= LUik * LUkj; // LU[i, k] LU[k, j] } } } } for (int i = 0; i < n; i++) { System.Diagnostics.Debug.Assert(LU.RowColIndexValues[i].ContainsKey(i)); } }
public static bool DoubleSolvePreconditionedBiCGSTAB( out double[] X, DoubleSparseMatrix A, double[] B, DoubleSparseMatrix LU, double convRatioTolerance) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; System.Diagnostics.Debug.Assert(B.Length == n); double convRatio = convRatioTolerance; double tolerance = convRatio; const int maxIter = IvyFEM.Linear.Constants.MaxIter; double[] r = new double[n]; double[] r0 = new double[n]; double[] x = new double[n]; double[] p = new double[n]; double[] s = new double[n]; double[] Mp = new double[n]; double[] Ms = new double[n]; int iter = 0; B.CopyTo(r, 0); double sqInvNorm0; { double sqNorm0 = IvyFEM.Lapack.Functions.ddot(r, r); if (sqNorm0 < IvyFEM.Constants.PrecisionLowerLimit) { convRatio = 0; X = x; System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); return(true); } sqInvNorm0 = 1.0 / sqNorm0; } r.CopyTo(r0, 0); r.CopyTo(p, 0); for (iter = 0; iter < maxIter; iter++) { DoubleSolveLU(out Mp, LU, p); double r0r = IvyFEM.Lapack.Functions.ddot(r0, r); double[] AMp = A * Mp; double alpha; { double denominator = IvyFEM.Lapack.Functions.ddot(r0, AMp); alpha = r0r / denominator; } s = IvyFEM.Lapack.Functions.daxpy(-alpha, AMp, r); DoubleSolveLU(out Ms, LU, s); double[] AMs = A * Ms; double omega; { double denominator = IvyFEM.Lapack.Functions.ddot(AMs, AMs); double numerator = IvyFEM.Lapack.Functions.ddot(s, AMs); omega = numerator / denominator; } x = IvyFEM.Lapack.Functions.daxpy(alpha, Mp, x); x = IvyFEM.Lapack.Functions.daxpy(omega, Ms, x); r = IvyFEM.Lapack.Functions.daxpy(-omega, AMs, s); { double sqNorm = IvyFEM.Lapack.Functions.ddot(r, r); if (sqNorm * sqInvNorm0 < tolerance * tolerance) { convRatio = Math.Sqrt(sqNorm * sqInvNorm0); X = x; System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); return(true); } } double beta; { double r0rPrev = r0r; r0r = IvyFEM.Lapack.Functions.ddot(r0, r); beta = (r0r * alpha) / (r0rPrev * omega); } p = IvyFEM.Lapack.Functions.daxpy(beta, p, r); p = IvyFEM.Lapack.Functions.daxpy(-beta * omega, AMp, p); } { double sqNormRes = IvyFEM.Lapack.Functions.ddot(r, r); convRatio = Math.Sqrt(sqNormRes * sqInvNorm0); System.Diagnostics.Debug.WriteLine("iter = " + iter + " norm: " + convRatio); X = x; } System.Diagnostics.Debug.WriteLine("Not converged"); return(false); }
public DoubleSparseMatrix(DoubleSparseMatrix src) { Copy(src); }
// Incomplete Cholesky Factorization public static DoubleSparseMatrix DoubleCalcIC(DoubleSparseMatrix A) { System.Diagnostics.Debug.Assert(A.RowLength == A.ColumnLength); int n = A.RowLength; DoubleSparseMatrix LU = new DoubleSparseMatrix(n, n); double[] D = new double[n]; D[0] = A[0, 0]; LU[0, 0] = 1.0; // L for (int i = 1; i < n; i++) { for (int j = 0; j <= (i - 1); j++) { if (!A.RowColIndexValues[i].ContainsKey(j)) { continue; } double tmp = A[i, j]; foreach (var pair in LU.RowColIndexValues[i]) { int k = pair.Key; double LUik = pair.Value; if (k >= 0 && k <= (j - 1)) { // } else { continue; } tmp -= LUik * LU[j, k] * D[k]; // LU[i, k] } LU[i, j] = (1.0 / D[j]) * tmp; } // i == jのとき { double tmp = A[i, i]; foreach (var pair in LU.RowColIndexValues[i]) { int k = pair.Key; double LUik = pair.Value; if (k >= 0 && k <= (i - 1)) { // } else { continue; } tmp -= LUik * LUik * D[k]; // LU[i, k] } D[i] = tmp; LU[i, i] = 1.0; } } // DL^T for (int i = 0; i < n; i++) { foreach (var pair in LU.RowColIndexValues[i]) { int j = pair.Key; double LUij = pair.Value; if (j >= 0 && j <= (i - 1)) { // } else { continue; } LU[j, i] = D[j] * LUij; // LU[i, j] } LU[i, i] = D[i]; } return(LU); }
public static bool OrderToDoubleBandMatrix( out DoubleSparseMatrix orderedA, out double[] orderedB, out int[] indexs, DoubleSparseMatrix A, double[] B) { orderedA = null; orderedB = null; indexs = null; // バンド幅を縮小する // 非0要素のパターンを取得 bool[,] matPattern = GetDoubleMatrixNonzeroPattern(A); int n = matPattern.GetLength(0); // subdiagonal、superdiagonalのサイズを取得する int iniSubdiaLength = 0; int iniSuperdiaLength = 0; { System.Diagnostics.Debug.WriteLine("/////initial Band Matrix info///////"); int rowcolLength; int subdiaLength; int superdiaLength; GetBandMatrixSubDiaSuperDia(matPattern, out rowcolLength, out subdiaLength, out superdiaLength); iniSubdiaLength = subdiaLength; iniSuperdiaLength = superdiaLength; } // 非0要素出現順に節点番号を格納 IList <int> newIndexs = new List <int>(); Queue <int> check = new Queue <int>(); int[] remain = new int[matPattern.GetLength(0)]; for (int i = 0; i < matPattern.GetLength(0); i++) { remain[i] = i; } while (newIndexs.Count < n) { for (int iRemain = 0; iRemain < remain.Length; iRemain++) { int i = remain[iRemain]; if (i == -1) { continue; } check.Enqueue(i); remain[iRemain] = -1; break; } while (check.Count > 0) { int i = check.Dequeue(); newIndexs.Add(i); for (int iRemain = 0; iRemain < remain.Length; iRemain++) { int j = remain[iRemain]; if (j == -1) { continue; } if (matPattern[i, j]) { check.Enqueue(j); remain[iRemain] = -1; } } } } System.Diagnostics.Debug.Assert(newIndexs.Count == n); DoubleSparseMatrix newA = new DoubleSparseMatrix(n, n); // 遅い //for (int row = 0; row < n; row++) //{ // for (int col = 0; col < n; col++) // { // newA[row, col] = A[newIndexs[row], newIndexs[col]]; // } //} int[] oldIndexs = new int[n]; for (int i = 0; i < n; i++) { oldIndexs[newIndexs[i]] = i; } for (int row = 0; row < n; row++) { foreach (var pair in A.RowColIndexValues[row]) { int col = pair.Key; double value = pair.Value; newA[oldIndexs[row], oldIndexs[col]] = value; } } // 改善できないこともあるのでチェックする bool improved = false; // 非0パターンを取得 bool[,] newMatPattern = GetDoubleMatrixNonzeroPattern(newA); // check { System.Diagnostics.Debug.WriteLine("///// orderd Band Matrix info ///////"); int rowcolLength; int subdiaLength; int superdiaLength; GetBandMatrixSubDiaSuperDia(newMatPattern, out rowcolLength, out subdiaLength, out superdiaLength); if (subdiaLength <= iniSubdiaLength && superdiaLength <= iniSuperdiaLength) { improved = true; } } if (improved) { // 置き換え orderedA = newA; indexs = newIndexs.ToArray(); orderedB = new double[n]; for (int row = 0; row < n; row++) { orderedB[row] = B[newIndexs[row]]; } } else { System.Diagnostics.Debug.WriteLine("band with not optimized!"); } return(improved); }