![]() |
Contains subroutines for sparse matrix assembly and the solver core. More...
Functions/Subroutines | |
subroutine, public | sparse |
Determines the position of jacobian entries in the cscf value array. More... | |
subroutine, public | netsolve (rhs, res) |
The solver core. More... | |
Variables | |
real(r_kind), dimension(:), allocatable, public | vals |
contains non-zero matrix elements (of the Jacobian) More... | |
integer, dimension(:), allocatable, public | rows |
rows(i) is the number of the row in the matrix that contains the i-th value in vals More... | |
integer, dimension(:), allocatable, public | pt_b |
pt_b(j) gives the index into vals that contains the FIRST non-zero element of column j of the matrix More... | |
integer, dimension(:), allocatable, public | pt_e |
pt_e(j)-pt_b(j) gives the index into vals that contains the LAST non-zero element of column j of the matrix More... | |
integer, dimension(:), allocatable, public | dia |
dia(j) gives the index into vals that contains the (j,j) diagonal element of the matrix More... | |
integer, dimension(:,:), allocatable, public | jind |
jind(i,j) contains the index of jacobian entry in vals More... | |
Contains subroutines for sparse matrix assembly and the solver core.
Edited:
subroutine, public pardiso_class::netsolve | ( | real(r_kind), dimension(net_size), intent(inout) | rhs, |
real(r_kind), dimension(net_size), intent(inout) | res | ||
) |
The solver core.
This subroutine calls the sparse PARDISO solver intel pardiso to solve equations in the following form:
\[ J\cdot \mathrm{res} = \mathrm{rhs} \]
with res being the unknown. In addition, this subroutine checks for zero entries in vals and adjusts the sparse format accordingly.
[in,out] | rhs | right-hand sides of the system |
[in,out] | res | resulting abundances |
Definition at line 246 of file pardiso_class.f90.
subroutine, public pardiso_class::sparse |
Determines the position of jacobian entries in the cscf value array.
This subroutine creates the basic sparse format of the jacobian. For this it writes a dummy matrix with entries = 1 at places in the matrix that are connected via nuclear reaction rates. Afterwards the matrix is stored in compressed sparse column format (cscf) into the arrays vals, rows, and pt_b. In cscf format, the matrix:
\[ X= \left( \begin{array}{cccc} 9 & 8 & 0 & 2\\ 0 & 0 & 0 & 5\\ 0 & 3 & 0 & 0\\ 5 & 0 & 4 & 9\end{array} \right). \]
will be stored like
\[ \mathrm{vals} =[9,5,8,3,4,2,5,9] \]
with vals storing the non-zero values of X,
\[ \mathrm{rows} = [1,4,1,3,4,1,2,4] \]
rows storing the position of a new row in vals, and
\[ \mathrm{pt\_b} = [1,3,5,6,9] \]
. the position of the columns of the non-zero entries.
Edited:
Definition at line 78 of file pardiso_class.f90.
integer, dimension(:), allocatable, public pardiso_class::dia |
dia(j) gives the index into vals that contains the (j,j) diagonal element of the matrix
Definition at line 34 of file pardiso_class.f90.
integer, dimension(:,:), allocatable, public pardiso_class::jind |
jind(i,j) contains the index of jacobian entry in vals
Definition at line 35 of file pardiso_class.f90.
integer, dimension(:), allocatable, public pardiso_class::pt_b |
pt_b(j) gives the index into vals that contains the FIRST non-zero element of column j of the matrix
Definition at line 32 of file pardiso_class.f90.
integer, dimension(:), allocatable, public pardiso_class::pt_e |
pt_e(j)-pt_b(j) gives the index into vals that contains the LAST non-zero element of column j of the matrix
Definition at line 33 of file pardiso_class.f90.
integer, dimension(:), allocatable, public pardiso_class::rows |
rows(i) is the number of the row in the matrix that contains the i-th value in vals
Definition at line 31 of file pardiso_class.f90.
real(r_kind), dimension(:), allocatable, public pardiso_class::vals |
contains non-zero matrix elements (of the Jacobian)
Definition at line 30 of file pardiso_class.f90.