SSIDS - Sparse Symmetric Indefinite Direct Solver
#include <spral_ssids.h> /* or <spral.h> for all packages */
Purpose
This package solves one or more sets of \(n\times n\) sparse symmetric equations \(AX = B\) using a multifrontal method. The following cases are covered:
1. \(A\) is indefinite. SSIDS computes the sparse factorization
where \(P\) is a permutation matrix, \(L\) is unit lower triangular, and \(D\) is block diagonal with blocks of size \(1 \times 1\) and \(2 \times 2\).
2. \(A\) is positive definite. SSIDS computes the sparse Cholesky factorization
where \(P\) is a permutation matrix and \(L\) is lower triangular. However, as SSIDS is designed primarily for indefinite systems, this may be slower than a dedicated Cholesky solver.
The code optionally supports hybrid computation using one or more NVIDIA GPUs.
SSIDS returns bit-compatible results.
An option exists to scale the matrix. In this case, the factorization of the scaled matrix \(\overline{A} = S A S\) is computed, where \(S\) is a diagonal scaling matrix.
Usage overview
Solving \(AX=B\) using SSIDS is a four stage process.
Call
spral_ssids_analyse()
to perform a symbolic factorization, stored in akeep.Call
spral_ssids_factor()
to perform a numeric factorization, stored in fkeep. More than one numeric factorization can refer to the same akeep.Call
spral_ssids_solve1()
orspral_ssids_solve()
to perform a solve with the factors. More than one solve can be performed with the same fkeep.Once all desired solutions have been performed, free memory with
spral_ssids_free()
.
In addition, advanced users may use the following routines:
spral_ssids_enquire_posdef()
andspral_ssids_enquire_indef()
return the diagonal entries of the factors and the pivot sequence.spral_ssids_alter()
allows altering the diagonal entries of the factors.
Note
Bit-compatibility: If used with bit-compatible BLAS and compiler options, this routine will return bit compatible results. That is, consecutive runs with the same data on the same machine produces exactly the same solution.
Basic Subroutines
Note
For the most efficient use of the package, CSC format should be used without checking.
-
void spral_ssids_default_options(struct spral_ssids_options *options)
Intialises members of options structure to default values.
- Parameters:
options – Structure to be initialised.
-
void spral_ssids_analyse(bool check, int n, int *order, const int64_t *ptr, const int *row, const double *val, void **akeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
Perform the analyse (symbolic) phase of the factorization for a matrix supplied in CSC format. The resulting symbolic factors stored in
akeep
should be passed unaltered in the following call tospral_ssids_factor()
.- Parameters:
check – if true, matrix data is checked. Out-of-range entries are dropped and duplicate entries are summed.
n – number of columns in \(A\).
order – may be NULL; otherwise must be an array of size n used on entry a user-supplied ordering (
options.ordering=0
). On return, the actual ordering used.ptr[n+1] – column pointers for \(A\) (see CSC format).
row[ptr[n]] – row indices for \(A\) (see CSC format).
val – may be NULL; otherwise must be an array of size ptr[n] containing non-zero values for \(A\) (see CSC format). Only used if a matching-based ordering is requested.
akeep – returns symbolic factorization, to be passed unchanged to subsequent routines.
options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).
Note
If a user-supplied ordering is used, it may be altered by this routine, with the altered version returned in order[]. This version will be equivalent to the original ordering, except that some supernodes may have been amalgamated, a topographic ordering may have been applied to the assembly tree and the order of columns within a supernode may have been adjusted to improve cache locality.
-
void spral_ssids_analyse_ptr32(bool check, int n, int *order, const int *ptr, const int *row, const double *val, void **akeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
As
spral_ssids_analyse()
, except ptr has typeint
. Provided for backwards comptability, users are encourage to use 64-bit ptr in new code.
-
void spral_ssids_analyse_coord(int n, int *order, int64_t ne, const int *row, const int *col, const double *val, void **akeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
As
spral_ssids_analyse()
, but for coordinate data. The variant parameters are:- Parameters:
ne – number of non-zero entries in \(A\).
row[ne] – row indices for \(A\) (see Coordinate format).
col[ne] – column indices for \(A\) (see Coordinate format).
-
void spral_ssids_factor(bool posdef, const int64_t *ptr, const int *row, const double *val, double *scale, void *akeep, void **fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
- Parameters:
posdef – true if matrix is positive-definite
ptr – may be NULL; otherwise a length n+1 array of column pointers for \(A\), only required if
akeep
was obtained by runningspral_ssids_analyse()
with check=true, in which case it must be unchanged since that call.row – may be NULL; otherwise a length ptr[n] array of row indices for \(A\), only required if
akeep
was obtained by runningspral_ssids_analyse()
with check=true, in which case it must be unchanged since that call.val[] – non-zero values for \(A\) in same format as for the call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.scale – may be NULL; otherwise a length n array for diagonal scaling. scale[i-1] contains entry \(S_ii\) of \(S\). Must be supplied by user on entry if
options.scaling=0
(user-supplied scaling). On exit, returns scaling used.akeep – symbolic factorization returned by preceding call to
ssids_analyse()
orssids_analyse_coord()
.fkeep – returns numeric factorization, to be passed unchanged to subsequent routines.
options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).
-
void spral_ssids_solve1(int job, double *x1, void *akeep, void *fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
Solve (for a single right-hand side) one of the following equations:
job
Equation solved
0
\(Ax=b\)
1
\(PLx=Sb\)
2
\(Dx=b\)
3
\((PL)^TS^{-1}x=b\)
4
\(D(PL)^TS^{-1}x=b\)
Recall \(A\) has been factorized as either:
\(SAS = (PL)(PL)^T~\) (positive-definite case); or
\(SAS = (PL)D(PL)^T\) (indefinite case).
- Parameters:
job – specifies equation to solve, as per above table.
x1[n] – right-hand side \(b\) on entry, solution \(x\) on exit.
akeep – symbolic factorization returned by preceding call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.fkeep – numeric factorization returned by preceding call to
spral_ssids_factor()
.options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).
-
void spral_ssids_solve(int job, int nrhs, double *x, int ldx, void *akeep, void *fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
Solve (for multiple right-hand sides) one of the following equations:
job
Equation solved
0
\(AX=B\)
1
\(PLX=SB\)
2
\(DX=B\)
3
\((PL)^TS^{-1}X=B\)
4
\(D(PL)^TS^{-1}X=B\)
Recall \(A\) has been factorized as either:
\(SAS = (PL)(PL)^T~\) (positive-definite case); or
\(SAS = (PL)D(PL)^T\) (indefinite case).
- Parameters:
job – specifies equation to solve, as per above table.
nrhs – number of right-hand sides.
x[ldx*nrhs] – right-hand sides \(B\) on entry, solutions \(X\) on exit. The i-th entry of right-hand side j is in position x[j*ldx+i].
ldx – leading dimension of x.
akeep – symbolic factorization returned by preceding call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.fkeep – numeric factorization returned by preceding call to
spral_ssids_factor()
.options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).
-
int spral_ssids_free_akeep(void **akeep)
Frees memory and resources associated with
akeep
.- Parameters:
akeep – symbolic factors to be freed.
- Returns:
0 on success, or a CUDA error code on failure.
-
int spral_ssids_free_fkeep(void **fkeep)
Frees memory and resources associated with
fkeep
.- Parameters:
fkeep – numeric factors to be freed.
- Returns:
0 on success, or a CUDA error code on failure.
-
int spral_ssids_free(void **akeep, void **fkeep)
Frees memory and resources associated with
akeep
andfkeep
.- Parameters:
akeep – symbolic factors to be freed.
fkeep – numeric factors to be freed.
- Returns:
0 on success, or a CUDA error code on failure.
Warning
The free routine(s) must be called by the user. Merely deallocating
akeep
or fkeep
, or allowing them to go out of scope will
result in memory leaks. akeep
should only be deallocated after all
associated numeric factorizations fkeep
have been freed.
Advanced subroutines
-
void spral_ssids_enquire_posdef(const void *akeep, const void *fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform, double *d)
Return the diagonal entries of the Cholesky factor.
- Parameters:
akeep – symbolic factorization returned by preceding call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.fkeep – numeric factorization returned by preceding call to
spral_ssids_factor()
.options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).d[n] – returns the diagonal of \(L\). d[i-1] stores the entry \(L_{ii}\).
-
void spral_ssids_enquire_indef(const void *akeep, const void *fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform, int *piv_order, double *d)
Return the pivot order and/or values of \(D\) of the Symmetric Indefinite Factorization.
- Parameters:
akeep – symbolic factorization returned by preceding call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.fkeep – numeric factorization returned by preceding call to
spral_ssids_factor()
.options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).piv_order – may be NULL; otherwise a length n array for the pivot order. On return, \(|\,\texttt{piv_order[i]}|\) gives the position of variable \(i\) in the pivot order.
d – may be NULL; otherwise a length 2*n array for the \(2\times2\) block diagonal of \(D\). d[2*(i-1)+0] stores \(D_{ii}\) and d[2*(i-1)+1] stores \(D_{(i+1)i}\).
-
void spral_ssids_alter(const double *d, const void *akeep, void *fkeep, const struct spral_ssids_options *options, struct spral_ssids_inform *inform)
Alter the entries of the diagonal factor \(D\) for a symmetric indefinite factorization. The pivot order remains the same.
- Parameters:
d[2*n] – New entries of \(D\). d[2*(i-1)+0] stores \(D_{ii}\) and d[2*(i-1)+1] stores \(D_{(i+1)i}\).
akeep – symbolic factorization returned by preceding call to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
.fkeep – numeric factorization returned by preceding call to
spral_ssids_factor()
.options – specifies algorithm options to be used (see
spral_ssids_options
).inform – returns information about the execution of the routine (see
spral_ssids_inform
).
Note: This routine is not compatabile with the option
options.presolve=1
.
Derived types
- struct spral_ssids_options
The derived data type ssids_options is used to specify the options used within
SSIDS
. The components, that are automatically given default values in the definition of the type, are:-
int array_base
Indexing base for arrays. Either 0 (C indexing) or 1 (Fortran indexing). Default is 0.
-
int print_level
Level of printing:
< 0
No printing.
= 0 (default)
Error and warning messages only.
= 1
As 0, plus basic diagnostic printing.
> 1
As 1, plus some additional diagnostic printing.
The default is 0.
-
int unit_diagnostics
Fortran unit number for diagnostics printing. Printing is suppressed if <0. The default is 6 (stdout).
-
int unit_error
Fortran unit number for printing of error messages. Printing is suppressed if <0. The default is 6 (stdout).
-
int unit_warning
Fortran unit number for printing of warning messages. Printing is suppressed if <0. The default is 6 (stdout).
-
int ordering
Ordering method to use in analyse phase:
0
User-supplied ordering is used (order argument to
spral_ssids_analyse()
orspral_ssids_analyse_coord()
).1 (default)
METIS ordering with default settings.
2
Matching-based elimination ordering is computed (the Hungarian algorithm is used to identify large off-diagonal entries. A restricted METIS ordering is then used that forces these on to the subdiagonal).
Note: This option should only be chosen for indefinite systems. A scaling is also computed that may be used in
spral_ssids_factor()
(seescaling
below).The default is 1.
-
int nemin
Supernode amalgamation threshold. Two neighbours in the elimination tree are merged if they both involve fewer than nemin eliminations. The default is used if nemin<1. The default is 32.
-
float gpu_perf_coeff
GPU performance coefficient. How many times faster a GPU is than CPU at factoring a subtree. Default is 1.0.
-
int scaling
Scaling algorithm to use:
<=0 (default)
No scaling (if
scale[]
is not present on call tospral_ssids_factor()
, or user-supplied scaling (ifscale[]
is present).=1
Compute using weighted bipartite matching via the Hungarian Algorithm (MC64 algorithm).
=2
Compute using a weighted bipartite matching via the Auction Algorithm (may be lower quality than that computed using the Hungarian Algorithm, but can be considerably faster).
=3
Use matching-based ordering generated during the analyse phase using
options.ordering=2
. The scaling will be the same as that generated withoptions.scaling=1
if the matrix values have not changed. This option will generate an error if a matching-based ordering was not used during analysis.>=4
Compute using the norm-equilibration algorithm of Ruiz (see SCALING - Sparse matrix scalings).
The default is 0.
-
int64_t small_subtree_threshold
Maximum number of flops in a subtree treated as a single task. See method section. The default is 4e6.
-
int cpu_block_size
Block size to use for parallelization of large nodes on CPU resources. Default is 256.
-
bool action
Continue factorization of singular matrix on discovery of zero pivot if true (a warning is issued), or abort if false. The default is true.
-
int pivot_method
Pivot method to be used on CPU, one of:
1
Aggressive a posteori pivoting. Cholesky-like communication pattern is used, but a single failed pivot requires restart of node factorization and potential recalculation of all uneliminated entries.
2 (default)
Block a posteori pivoting. A failed pivot only requires recalculation of entries within its own block column.
3
Threshold partial pivoting. Not parallel.
Default is 2.
-
double small
Threshold below which an entry is treated as equivalent to 0.0. The default is 1e-20.
-
double u
Relative pivot threshold used in symmetric indefinite case. Values outside of the range \([0,0.5]\) are treated as the closest value in that range. The default is 0.01.
-
int array_base
- struct spral_ssids_inform
Used to return information about the progress and needs of the algorithm.
-
int64_t cpu_flops
Number of flops performed on CPU
-
int cublas_error
CUBLAS error code in the event of a CUBLAS error (0 otherwise).
-
int cuda_error
CUDA error code in the event of a CUDA error (0 otherwise). Note that due to asynchronous execution, CUDA errors may not be reported by the call that caused them.
-
int flag
Exit status of the algorithm (see table below).
-
int64_t gpu_flops
Number of flops performed on GPU
-
int matrix_dup
Number of duplicate entries encountered (if
spral_ssids_analyse()
called with check=true, or any call tospral_ssids_analyse_coord()
).
-
int matrix_missing_diag
Number of diagonal entries without an explicit value (if
spral_ssids_analyse()
called with check=true, or any call tospral_ssids_analyse_coord()
).
-
int matrix_outrange
Number of out-of-range entries encountered (if
spral_ssids_analyse()
called with check=true, or any call tospral_ssids_analyse_coord()
).
-
int matrix_rank
(Estimated) rank (structural after analyse phase, numerical after factorize phase).
-
int maxdepth
Maximum depth of the assembly tree.
-
int maxfront
Maximum front size (without pivoting after analyse phase, with pivoting after factorize phase).
-
int maxsupernode
Maximum supernode size (without pivoting after analyse phase, with pivoting after factorize phase).
-
int num_delay
Number of delayed pivots. That is, the total number of fully-summed variables that were passed to the father node because of stability considerations. If a variable is passed further up the tree, it will be counted again.
-
int64_t num_factor
Number of entries in \(L\) (without pivoting after analyse phase, with pivoting after factorize phase).
-
int64_t num_flops
Number of floating-point operations for Cholesky factorization (indefinte needs slightly more). Without pivoting after analyse phase, with pivoting after factorize phase.
-
int num_neg
Number of negative eigenvalues of the matrix \(D\) after factorize phase.
-
int num_sup
Number of supernodes in assembly tree.
-
int num_two
Number of \(2 \times 2\) pivots used by the factorization (i.e. in the matrix \(D\)).
-
int stat
Fortran allocation status parameter in event of allocation error (0 otherwise).
inform.flag
Return status
0
Success.
-1
Error in sequence of calls (may be caused by failure of a preceding call).
-2
n<0 or ne<1.
-3
Error in ptr[].
-4
CSC format: All variable indices in one or more columns are out-of-range.
Coordinate format: All entries are out-of-range.
-5
Matrix is singular and options.action=false
-6
Matrix found not to be positive definite.
-7
ptr[] and/or row[] not present, but required as
spral_ssids_analyse()
was called with check=false.-8
options.ordering out of range, or options.ordering=0 and order parameter not provided or not a valid permutation.
-9
options.ordering=-2 but val[] was not supplied.
-10
ldx<n or nrhs<1.
-11
job is out-of-range.
-13
Called
spral_ssids_enquire_posdef()
on indefinite factorization.-14
Called
spral_ssids_enquire_indef()
on positive-definite factorization.-15
options.scaling=3 but a matching-based ordering was not performed during analyse phase.
-50
Allocation error. If available, the stat parameter is returned in inform.stat.
-51
CUDA error. The CUDA error return value is returned in inform.cuda_error.
-52
CUBLAS error. The CUBLAS error return value is returned in inform.cublas_error.
-53
OpenMP cancellation is disabled. Please set the environment variable OMP_CANCELLATION=true.
+1
Out-of-range variable indices found and ignored in input data. inform.matrix_outrange is set to the number of such entries.
+2
Duplicate entries found and summed in input data. inform.matrix_dup is set to the number of such entries.
+3
Combination of +1 and +2.
+4
One or more diagonal entries of \(A\) are missing.
+5
Combination of +4 and +1 or +2.
+6
Matrix is found be (structurally) singular during analyse phase. This will overwrite any of the above warning flags.
+7
Matrix is found to be singular during factorize phase.
+8
Matching-based scaling found as side-effect of matching-based ordering ignored (consider setting options.scaling=3).
+50
OpenMP processor binding is disabled. Consider setting the environment variable OMP_PROC_BIND=true (this may affect performance on NUMA systems).
-
int64_t cpu_flops
Example
Suppose we wish to factorize the matrix
and then solve for the right-hand side
The following code may be used.
/* examples/C/ssids.c - Example code for SPRAL_SSIDS package */
#include "spral.h"
#include <stdlib.h>
#include <stdint.h>
#include <stdio.h>
int main(void) {
/* Derived types */
void *akeep, *fkeep;
struct spral_ssids_options options;
struct spral_ssids_inform inform;
/* Initialize derived types */
akeep = NULL; fkeep = NULL; /* Important that these are NULL to start with */
spral_ssids_default_options(&options);
options.array_base = 1; /* Need to set to 1 if using Fortran 1-based indexing */
/* Data for matrix:
* ( 2 1 )
* ( 1 4 1 1 )
* ( 1 3 2 )
* ( 2 -1 )
* ( 1 2 ) */
bool posdef = false;
int n = 5;
int64_t ptr[] = { 1, 3, 6, 8, 9, 10 };
int row[] = { 1, 2, 2, 3, 5, 3, 4, 4, 5 };
double val[] = { 2.0, 1.0, 4.0, 1.0, 1.0, 3.0, 2.0, -1.0, 2.0 };
/* The right-hand side with solution (1.0, 2.0, 3.0, 4.0, 5.0) */
double x[] = { 4.0, 17.0, 19.0, 2.0, 12.0 };
/* Perform analyse and factorise with data checking */
bool check = true;
spral_ssids_analyse(check, n, NULL, ptr, row, NULL, &akeep, &options,
&inform);
if(inform.flag<0) {
spral_ssids_free(&akeep, &fkeep);
exit(1);
}
spral_ssids_factor(posdef, NULL, NULL, val, NULL, akeep, &fkeep, &options,
&inform);
if(inform.flag<0) {
spral_ssids_free(&akeep, &fkeep);
exit(1);
}
/* Solve */
spral_ssids_solve1(0, x, akeep, fkeep, &options, &inform);
if(inform.flag<0) {
spral_ssids_free(&akeep, &fkeep);
exit(1);
}
printf("The computed solution is:\n");
for(int i=0; i<n; i++) printf(" %18.10e", x[i]);
printf("\n");
/* Determine and print the pivot order */
int piv_order[5];
spral_ssids_enquire_indef(akeep, fkeep, &options, &inform, piv_order, NULL);
printf("Pivot order:");
for(int i=0; i<n; i++) printf(" %5d", piv_order[i]);
printf("\n");
int cuda_error = spral_ssids_free(&akeep, &fkeep);
if(cuda_error!=0) exit(1);
return 0;
}
This produces the following output:
The computed solution is:
1.0000000000e+00 2.0000000000e+00 3.0000000000e+00 4.0000000000e+00 5.0000000000e+00
Pivot order: 2 3 0 -1 1
Driver Program
SSIDS ships with a driver program spral_ssids
that allows reading a
matrix in Rutherford-Boeing format specified as a command-line argument and
factorizing it. There are a number of other command-line arguments that
configure the factorization.
- spral_ssids
SSIDS driver program.
- Parameters:
filename –
Rutherford-Boeing matrix filename (default if not specified is matrix.rb).
- --scale=none
use no scaling (the default).
- --scale=mc64
use the Hungarian scaling algorithm (as in MC64).
- --scale=auction
use the Auction scaling algorithm.
- --scale=mc77
use the norm-equilibration scaling algorithm (as in MC77).
- --ordering=mc64-metis
use matching-based ordering and scaling (scale is overwritten).
- --force-posdef
force the matrix to be positive definite
- --posdef
assume the matrix is positive definite.
- --time-scaling
time the scaling routine.
- --nrhs
set the number of right-hand sides [integer,default=1].
- --nemin
set the supernode amalgamation threshold [integer,default=32].
- --u
set the relative pivot threshold used in the symmetric indefinite case [real,default=0.01].
- --max-load-inbalance
set the maximum permissible load inbalance for leaf subtree allocations [real,default=1.2].
- --pivot-method=app-aggressive
use aggressive a posteori pivoting.
- --pivot-method=app-block
use block a posteori pivoting (the default).
- --pivot-method=tpp
use threshold partial pivoting.
- --flat-topology
force a flat machine topology (the default).
- --no-flat-topology
use the actual machine topology.
- --disable-gpu
don’t use an NVIDIA GPU if present.
- --min-gpu-work
set the minimum number of flops in a subtree before scheduling on GPU [integer(long),default=5e9].
- --gpu-perf-coeff
set the GPU performance coefficient (how many times faster a GPU is than CPU at factoring a subtree) [real,default=1.0].
- --small-subtree-threshold
set the maximum number of flops in a subtree treated as a single task [integer(long),default=4e6].
- --cpu-block-size
set the block size to use for parallelization of large nodes on CPU resources [integer ,default=256].
- --no-ignore-numa
don’t treat all CPUs and GPUs as belonging to a single NUMA region (which is the default).
- --ngpus
set the number of NVIDIA GPUs to use [integer,default=0].
For example, to use auction scaling with two right-hand sides on the linverse.rb matrix:
./spral_ssids linverse.rb --scale=auction --nrhs 2
This produces output similar to the following:
The computed solution is:
Set scaling to Auction
solving for 2 right-hand sides
Reading 'linverse.rb'...
ok
Number of CUDA devices: 0
Forcing topology to 32
Using 0 GPUs
Used order 1
ok
Analyse took 5.20000011E-02
Predict nfact = 3.03E+05
Predict nflop = 9.25E+06
nparts 1
cpu_fl 9.25E+06
gpu_fl 0.00E+00
Factorize...
ok
Factor took 1.20000001E-02
Solve...
ok
Solve took 1.00000005E-03
number bad cmp = 0
fwd error || ||_inf = 3.8014036363165360E-013
bwd error scaled = 4.2549737582555113E-015 4.2549737582555113E-015
cmp: SMFCT
anal: 0.05
fact: 0.01
afact: 3.03E+05
aflop: 9.25E+06
nfact: 3.03E+05
nflop: 9.25E+06
delay: 0
inerti 2838 0 9161
2x2piv 5502
maxfro 55
maxsup 52
not_fi 0
not_se 0
Method
Partition of work across available resources
Once the ordering has been determined and the assembly tree determined in the analyse phase, the tree is broken into a number of leaf subtrees rooted at a single node, leaving a root subtree/forest consisting of all remaining nodes above those. Each leaf subtree is pre-assigned to a particular NUMA region or GPU for the factorization phase. Details of the algorithm used for finding these subtrees and their assignment can be found in the paper [2].
The factorization phase has two steps. In the first, leaf subtrees are factorized in parallel on their assigned resources. Once all leaf subtrees are factored, the second phase begins where all CPU resources cooperate to factorize the root subtree.
At present the solve phase is performed in serial.
Data checking
If check is set to .true. on the call to spral_ssids_analyse()
or if
spral_ssids_analyse_coord()
is called, the user-supplied matrix data is
checked for errors. The cleaned integer matrix data (duplicates are
summed and out-of-range indices discarded) is stored in akeep. The use
of checking is optional on a call to spral_ssids_analyse()
as it incurs both
time and memory overheads. However, it is recommended since the
behaviour of the other routines in the package is unpredictable if
duplicates and/or out-of-range variable indices are entered.
If the user has supplied an elimination order it is checked for errors.
Otherwise, an elimination order is generated by the package. The
elimination order is used to construct an assembly tree. On exit from
spral_ssids_analyse()
(and spral_ssids_analyse_coord()
), order[] is
set so that order[i] holds the position of variable \(i\) in the
elimination order. If an ordering was supplied by the user, this order may
differ, but will be equivalent in terms of fill-in.
Factorization performed
The factorization performed depends on the value of posdef.
If posdef=true, a Cholesky factorization is performed:
Pivoting is not performed, so \(P\) is the permutation determined in the analysis phase. \(S\) is a diagonal scaling matrix.
If posdef=false, a symmetric indefinite factorization is performed:
Pivoting is performed, so \(P\) may differ from the permutation determined in the analysis phase, though it is kept as close as possible to minimize fill. The exact pivoting algorithm varies depending on the particular kernel the algorithm chooses to employ.
Full details of the algorithms used are provided in [1] for GPUs and [2] for CPUs.
Small Leaf Subtrees
For subtrees allocated to run on the CPU, the factorization of small nodes near
the leaves of the tree can be amalgamated into a single parallel task (normally
each would be treated as its own OpenMP task to be scheduled). This can reduce
scheduling overheads, especially on small problems. If the total number of
operations for a subtree root at a given node is less than
options.small_subtree_threshold
,
that subtree is treated as a single task.