regrid_solvers module reference¶
Solvers of linear systems.
Functions/Subroutines¶
Solve the linear system AX = R by Gaussian elimination. |
|
Solve the linear system AX = R by Gaussian elimination. |
|
Solve the tridiagonal system AX = R. |
|
Solve the tridiagonal system AX = R. |
Detailed Description¶
Date of creation: 2008.06.12 L. White.
This module contains solvers of linear systems. These routines have now been updated for greater efficiency, especially in special cases.
Function/Subroutine Documentation¶
-
subroutine
regrid_solvers/
solve_linear_system
(A, R, X, N, answer_date)¶ Solve the linear system AX = R by Gaussian elimination.
This routine uses Gauss’s algorithm to transform the system’s original matrix into an upper triangular matrix. Back substitution yields the answer. The matrix A must be square, with the first index varing down the column.
- Parameters:
n :: [in] The size of the system
a :: [inout] The matrix being inverted in arbitrary units [A] on input, but internally modified to become nondimensional during the solver.
r :: [inout] system right-hand side in arbitrary units [A B] on input, but internally modified to have units of [B] during the solver
x :: [inout] solution vector in arbitrary units [B]
answer_date :: [in] The vintage of the expressions to use
- Call to:
-
subroutine
regrid_solvers/
linear_solver
(N, A, R, X)¶ Solve the linear system AX = R by Gaussian elimination.
This routine uses Gauss’s algorithm to transform the system’s original matrix into an upper triangular matrix. Back substitution then yields the answer. The matrix A must be square, with the first index varing along the row.
- Parameters:
n :: [in] The size of the system
a :: [inout] The matrix being inverted in arbitrary units [A] on input, but internally modified to become nondimensional during the solver.
r :: [inout] system right-hand side in [A B] on input, but internally modified to have units of [B] during the solver
x :: [inout] solution vector [B]
- Call to:
-
subroutine
regrid_solvers/
solve_tridiagonal_system
(Al, Ad, Au, R, X, N, answer_date)¶ Solve the tridiagonal system AX = R.
This routine uses Thomas’s algorithm to solve the tridiagonal system AX = R. (A is made up of lower, middle and upper diagonals)
- Parameters:
n :: [in] The size of the system
ad :: [in] Matrix center diagonal in arbitrary units [A]
al :: [in] Matrix lower diagonal [A]
au :: [in] Matrix upper diagonal [A]
r :: [in] system right-hand side in arbitrary units [A B]
x :: [out] solution vector in arbitrary units [B]
answer_date :: [in] The vintage of the expressions to use
- Called from:
regrid_edge_values::edge_slopes_implicit_h3
regrid_edge_values::edge_slopes_implicit_h5
regrid_edge_values::edge_values_implicit_h4
regrid_edge_values::edge_values_implicit_h6
-
subroutine
regrid_solvers/
solve_diag_dominant_tridiag
(Al, Ac, Au, R, X, N)¶ Solve the tridiagonal system AX = R.
This routine uses a variant of Thomas’s algorithm to solve the tridiagonal system AX = R, in a form that is guaranteed to avoid dividing by a zero pivot. The matrix A is made up of lower (Al) and upper diagonals (Au) and a central diagonal Ad = Ac+Al+Au, where Al, Au, and Ac are all positive (or negative) definite. However when Ac is smaller than roundoff compared with (Al+Au), the answers are prone to inaccuracy.
- Parameters:
n :: [in] The size of the system
ac :: [in] Matrix center diagonal offset from Al + Au in arbitrary units [A]
al :: [in] Matrix lower diagonal [A]
au :: [in] Matrix upper diagonal [A]
r :: [in] system right-hand side in arbitrary units [A B]
x :: [out] solution vector in arbitrary units [B]
- Called from:
regrid_edge_values::edge_slopes_implicit_h3
regrid_edge_values::edge_values_implicit_h4