FE-Project
|
Module common / Linear algebra. More...
Data Types | |
interface | linalgebra_solvelineq |
Solve a linear equation Ax=b using a direct solver. More... |
Functions/Subroutines | |
real(rp) function, dimension(size(a, 1), size(a, 2)), public | linalgebra_inv (a) |
Calculate a inversion of matrix A. | |
subroutine, public | linalgebra_lu (a_lu, ipiv) |
Perform LU factorization. | |
subroutine, public | linalgebra_solvelineq_bndmat (a, b, ipiv, n, kl, ku, nrhs, use_lapack) |
Calculate the solution of linear equations with a band matrix. | |
subroutine, public | linalgebra_solvelineq_gmres (a, b, x, m, restart_num, conv_crit) |
Solve a linear equation Ax=b using the Generalized Minimal Residual solver (GMRES) |
Module common / Linear algebra.
real(rp) function, dimension(size(a,1),size(a,2)), public scale_linalgebra::linalgebra_inv | ( | real(rp), dimension(:,:), intent(in) | a | ) |
Calculate a inversion of matrix A.
[in] | a | Matrix A |
Definition at line 76 of file scale_linalgebra.F90.
Referenced by scale_element_base::elementbase_construct_massmat(), get_pmatd_lu(), scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().
subroutine, public scale_linalgebra::linalgebra_lu | ( | real(rp), dimension(:,:), intent(inout) | a_lu, |
integer, dimension(size(a_lu,1)), intent(out) | ipiv ) |
Perform LU factorization.
Definition at line 173 of file scale_linalgebra.F90.
Referenced by scale_atm_dyn_dgm_nonhydro3d_hevi_gmres::atm_dyn_dgm_nonhydro3d_hevi_cal_vi(), get_pmatd_lu(), and scale_atm_dyn_dgm_hydrostatic::hydrostaic_build_rho_xyz::hydrostaic_build_rho_xyz_moist().
subroutine, public scale_linalgebra::linalgebra_solvelineq_bndmat | ( | real(rp), dimension(2*kl+ku+1,n), intent(inout) | a, |
real(rp), dimension(n,nrhs), intent(inout) | b, | ||
integer, dimension(n), intent(out) | ipiv, | ||
integer, intent(in) | n, | ||
integer, intent(in) | kl, | ||
integer, intent(in) | ku, | ||
integer, intent(in) | nrhs, | ||
logical, intent(in), optional | use_lapack ) |
Calculate the solution of linear equations with a band matrix.
In this module, we treat a linear equations system, A * X = B where A is a band matrix of order N with KL subdiagonals and KU superdiagonals, and X and B are matrices whose size is N * NRHS.
A | Matrix A in band storage same as in dgbsv of LAPACK |
b | Right hand side matrix B with N-by-RHS |
ipiv | Pivot indices that define the permutation matrix |
KL | Number of subdiagonals within the band of A |
KU | Number of superdiagonals within the band of A |
NRHS | Number of right hand sides |
use_lapack | Flag whether LAPACK is used |
Definition at line 208 of file scale_linalgebra.F90.
Referenced by scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_vi_asis(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi_asis(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform::atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform_cal_vi().
subroutine, public scale_linalgebra::linalgebra_solvelineq_gmres | ( | type(sparsemat), intent(in) | a, |
real(rp), dimension(:), intent(in) | b, | ||
real(rp), dimension(:), intent(inout) | x, | ||
integer, intent(in) | m, | ||
integer, intent(in) | restart_num, | ||
real(rp), intent(in) | conv_crit ) |
Solve a linear equation Ax=b using the Generalized Minimal Residual solver (GMRES)
[in] | a | Object managing sparse matrix A |
[in] | b | Right hand side vector |
[in,out] | x | Solution vector |
[in] | restart_num | Number of times by which restarting is performed if convergence condition is not satisfied |
[in] | conv_crit | Threshold for finishing the iteration solver |
Definition at line 320 of file scale_linalgebra.F90.
References scale_sparsemat::sparsemat_storage_typeid_csr.