SpM Handbook 1.2.4
|
#include "common.h"
#include <lapacke.h>
#include <cblas.h>
#include "frobeniusupdate.h"
Go to the source code of this file.
#define | LAPACKE_zlassq_work(_n_, _X_, _incx_, _scale_, _sumsq_) __spm_zlassq( (_n_), (_X_), (_incx_), (_scale_), (_sumsq_) ) |
TODO. | |
static int | __spm_zlassq (spm_int_t n, spm_complex64_t *x, spm_int_t incx, double *scale, double *sumsq) |
Updates the values scale and sumsq such that. | |
static void | z_spm_frobenius_elt_sym_diag (spm_int_t dofs, const double *valptr, double *data) |
Compute the Frobenius norm of a diagonal element within a symmetric/hermitian matrix with column/row major storage. | |
static void | z_spm_frobenius_elt_sym_offd (spm_int_t nbelts, const double *valptr, double *data) |
Compute the Frobenius norm of an off-diagonal element matrix in the symmetric/hermitian case. | |
static void | z_spm_frobenius_elt_sym (spm_int_t row, spm_int_t dofi, spm_int_t col, spm_int_t dofj, const spm_complex64_t *valptr, double *data) |
Compute the Frobenius norm of any element matrix in the symmetric/hermitian case. | |
static void | z_spmFrobeniusNorm_csc (const spmatrix_t *spm, double *data) |
Compute the Frobenius norm of a symmetrix/hermitian CSC matrix. | |
static void | z_spmFrobeniusNorm_csr (const spmatrix_t *spm, double *data) |
Compute the Frobenius norm of a symmetrix/hermitian CSR matrix. | |
static void | z_spmFrobeniusNorm_ijv (const spmatrix_t *spm, double *data) |
Compute the Frobenius norm of a symmetrix/hermitian IJV matrix. | |
double | z_spmFrobeniusNorm (const spmatrix_t *spm) |
Compute the Frobenius norm of the given spm structure. | |
double | z_spmMaxNorm (const spmatrix_t *spm) |
Compute the Max norm of the given spm structure. | |
static void | z_spm_oneinf_elt_sym_diag (spm_int_t row, spm_int_t dofi, const spm_complex64_t *valptr, double *sumtab) |
Compute the sum array for the one/inf norms of a diagonal element within a symmetric/hermitian matrix with column/row major storage. | |
static void | z_spm_oneinf_elt_gen_A (spm_int_t dofi, spm_int_t dofj, const spm_complex64_t *valptr, double *sumtab) |
Compute the sum array for the one/inf norms of a general element. | |
static void | z_spm_oneinf_elt_gen_B (spm_int_t dofi, spm_int_t dofj, const spm_complex64_t *valptr, double *sumtab) |
Compute the sum array for the one/inf norms of a general element. | |
static void | z_spm_oneinf_elt_gen_AB (spm_int_t row, spm_int_t dofi, spm_int_t col, spm_int_t dofj, const spm_complex64_t *valptr, double *sumtab) |
Compute the sum array for both the one and inf norms for the off-diagonal elements of symmetric/hermitian element matrices. | |
static void | z_spm_oneinf_elt_gen (spm_layout_t layout, spm_int_t row, spm_int_t dofi, spm_int_t col, spm_int_t dofj, const spm_complex64_t *valptr, spm_normtype_t ntype, double *sumtab) |
Compute the sum array for the one and inf norms of a general element. | |
static void | z_spm_oneinf_elt_sym_offd (spm_layout_t layout, spm_int_t row, spm_int_t dofi, spm_int_t col, spm_int_t dofj, const spm_complex64_t *valptr, double *sumtab) |
Compute the sum array for both the one and inf norms for the off-diagonal elements of symmetric/hermitian element matrices in either column or row major layout. | |
static void | z_spm_oneinf_elt (spm_mtxtype_t mtxtype, spm_layout_t layout, spm_int_t row, spm_int_t dofi, spm_int_t col, spm_int_t dofj, const spm_complex64_t *valptr, spm_normtype_t ntype, double *sumtab) |
Compute the sum array for the one/inf norm for an element matrix. | |
static void | z_spmOneInfNorm_csc (spm_normtype_t ntype, const spmatrix_t *spm, double *sumtab) |
Compute the one/inf norm of an spm CSC structure. | |
static void | z_spmOneInfNorm_csr (spm_normtype_t ntype, const spmatrix_t *spm, double *sumtab) |
Compute the one/inf norm of an spm CSR structure. | |
static void | z_spmOneInfNorm_ijv (spm_normtype_t ntype, const spmatrix_t *spm, double *sumtab) |
Compute the one/inf norm of an spm IJV structure. | |
static double | z_spmOneInfNorm (spm_normtype_t ntype, const spmatrix_t *spm) |
Compute the one/inf norm of the given spm structure given by the maximum row sum. | |
double | z_spmNorm (spm_normtype_t ntype, const spmatrix_t *spm) |
Compute the norm of an spm matrix. | |
double | z_spmNormMat (spm_normtype_t ntype, const spmatrix_t *spm, spm_int_t n, const spm_complex64_t *A, spm_int_t lda) |
Compute the norm of a dense matrix that follows the distribution of an spm matrix. | |
SParse Matrix package norm routine.
Definition in file z_spm_norm.c.
#define LAPACKE_zlassq_work | ( | _n_, | |
_X_, | |||
_incx_, | |||
_scale_, | |||
_sumsq_ | |||
) | __spm_zlassq( (_n_), (_X_), (_incx_), (_scale_), (_sumsq_) ) |
TODO.
Definition at line 88 of file z_spm_norm.c.
|
static |
Updates the values scale and sumsq such that.
( scale**2 )*sumsq = x( 1 )**2 +...+ x( n )**2 + ( scale**2 )*sumsq,
This routine is inspired from LAPACK zlassq function.
[in] | n | The number of elements in the vector |
[in] | x | The vector of size abs(n * incx) |
[in] | incx | The increment between two elements in the vector x. |
[in,out] | scale | On entry, the former scale On exit, the update scale to take into account the value |
[in,out] | sumsq | On entry, the former sumsq On exit, the update sumsq to take into account the value |
Definition at line 63 of file z_spm_norm.c.
|
static |
Compute the Frobenius norm of a diagonal element within a symmetric/hermitian matrix with column/row major storage.
Note that column major is using the low triangular part only of the diagonal element matrices, and row major, by symmetry, is using only the upper triangular part.
The comments in the code are made for column major storage.
[in] | dofs | TODO |
[in] | valptr | TODO |
[in,out] | data | TODO |
Definition at line 152 of file z_spm_norm.c.
Referenced by z_spm_frobenius_elt_sym().
|
static |
Compute the Frobenius norm of an off-diagonal element matrix in the symmetric/hermitian case.
[in] | nbelts | TODO |
[in] | valptr | TODO |
[in,out] | data |
Definition at line 206 of file z_spm_norm.c.
Referenced by z_spm_frobenius_elt_sym().
|
static |
Compute the Frobenius norm of any element matrix in the symmetric/hermitian case.
[in] | row | TODO |
[in] | dofi | TODO |
[in] | col | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in,out] | data | TODO |
Definition at line 252 of file z_spm_norm.c.
References z_spm_frobenius_elt_sym_diag(), and z_spm_frobenius_elt_sym_offd().
Referenced by z_spmFrobeniusNorm_csc(), z_spmFrobeniusNorm_csr(), and z_spmFrobeniusNorm_ijv().
|
static |
Compute the Frobenius norm of a symmetrix/hermitian CSC matrix.
[in] | spm | The spm from which the norm need to be computed. |
[in,out] | data | TODO |
Definition at line 283 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::replicated, spmatrix_s::rowptr, SpmComplex64, SpmCSC, spmatrix_s::values, and z_spm_frobenius_elt_sym().
Referenced by z_spmFrobeniusNorm().
|
static |
Compute the Frobenius norm of a symmetrix/hermitian CSR matrix.
[in] | spm | The spm from which the norm need to be computed. |
[in,out] | data | TODO |
Definition at line 350 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::loc2glob, spmatrix_s::n, spmatrix_s::replicated, spmatrix_s::rowptr, SpmComplex64, SpmCSR, spmatrix_s::values, and z_spm_frobenius_elt_sym().
Referenced by z_spmFrobeniusNorm().
|
static |
Compute the Frobenius norm of a symmetrix/hermitian IJV matrix.
[in] | spm | The spm from which the norm need to be computed. |
[in,out] | data | TODO |
Definition at line 418 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::flttype, spmatrix_s::fmttype, spmatrix_s::nnz, spmatrix_s::rowptr, SpmComplex64, SpmIJV, spmatrix_s::values, and z_spm_frobenius_elt_sym().
Referenced by z_spmFrobeniusNorm().
double z_spmFrobeniusNorm | ( | const spmatrix_t * | spm | ) |
Compute the Frobenius norm of the given spm structure.
||A|| = sqrt( sum( a_ij ^ 2 ) )
[in] | spm | The spm from which the norm need to be computed. |
Definition at line 480 of file z_spm_norm.c.
References spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::fmttype, spmatrix_s::mtxtype, spmatrix_s::nnzexp, spmatrix_s::replicated, SpmCSC, SpmCSR, SpmGeneral, SpmIJV, spmatrix_s::values, z_spmFrobeniusNorm_csc(), z_spmFrobeniusNorm_csr(), and z_spmFrobeniusNorm_ijv().
Referenced by z_spmNorm().
double z_spmMaxNorm | ( | const spmatrix_t * | spm | ) |
Compute the Max norm of the given spm structure.
||A|| = max( abs(a_ij) )
[in] | spm | The spm from which the norm need to be computed. |
Definition at line 543 of file z_spm_norm.c.
References spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::nnzexp, spmatrix_s::replicated, and spmatrix_s::values.
Referenced by z_spmNorm().
|
static |
Compute the sum array for the one/inf norms of a diagonal element within a symmetric/hermitian matrix with column/row major storage.
Note that column major is using the low triangular part only of the diagonal element matrices, and row major, by symmetry, is using only the upper triangular part.
The comments in the code are made for column major storage.
[in] | row | TODO |
[in] | dofi | TODO |
[in] | valptr | TODO |
[in,out] | sumtab | TODO |
Definition at line 591 of file z_spm_norm.c.
Referenced by z_spm_oneinf_elt().
|
static |
Compute the sum array for the one/inf norms of a general element.
We can observe two cases A and B;
| | | | | | |________| | | | | | | |A| OR |________| |B| |__|__|__| |_| |________| |_|
|___B____| |___A____|
| One Norm | Inf norm |
----------—+-------—+-------—+ Column Major | B | A | ----------—+-------—+-------—+ Row Major | A | B | ----------—+-------—+-------—+
[in] | dofi | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in,out] | sumtab | TODO |
Definition at line 659 of file z_spm_norm.c.
Referenced by z_spm_oneinf_elt_gen().
|
static |
Compute the sum array for the one/inf norms of a general element.
[in] | dofi | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in,out] | sumtab | TODO |
Definition at line 698 of file z_spm_norm.c.
Referenced by z_spm_oneinf_elt_gen().
|
static |
Compute the sum array for both the one and inf norms for the off-diagonal elements of symmetric/hermitian element matrices.
[in] | row | TODO |
[in] | dofi | TODO |
[in] | col | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in,out] | sumtab | TODO |
Definition at line 744 of file z_spm_norm.c.
Referenced by z_spm_oneinf_elt_sym_offd().
|
static |
Compute the sum array for the one and inf norms of a general element.
[in] | layout | TODO |
[in] | row | TODO |
[in] | dofi | TODO |
[in] | col | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in] | ntype | TODO |
[in,out] | sumtab | TODO |
Definition at line 799 of file z_spm_norm.c.
References SpmColMajor, SpmInfNorm, SpmOneNorm, z_spm_oneinf_elt_gen_A(), and z_spm_oneinf_elt_gen_B().
Referenced by z_spm_oneinf_elt().
|
static |
Compute the sum array for both the one and inf norms for the off-diagonal elements of symmetric/hermitian element matrices in either column or row major layout.
[in] | layout | TODO |
[in] | row | TODO |
[in] | dofi | TODO |
[in] | col | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in,out] | sumtab | TODO |
Definition at line 860 of file z_spm_norm.c.
References SpmColMajor, and z_spm_oneinf_elt_gen_AB().
Referenced by z_spm_oneinf_elt().
|
static |
Compute the sum array for the one/inf norm for an element matrix.
[in] | mtxtype | TODO |
[in] | layout | TODO |
[in] | row | TODO |
[in] | dofi | TODO |
[in] | col | TODO |
[in] | dofj | TODO |
[in] | valptr | TODO |
[in] | ntype | TODO |
[in,out] | sumtab | TODO |
Definition at line 912 of file z_spm_norm.c.
References SpmGeneral, z_spm_oneinf_elt_gen(), z_spm_oneinf_elt_sym_diag(), and z_spm_oneinf_elt_sym_offd().
Referenced by z_spmOneInfNorm_csc(), z_spmOneInfNorm_csr(), and z_spmOneInfNorm_ijv().
|
static |
Compute the one/inf norm of an spm CSC structure.
[in] | ntype | TODO |
[in] | spm | TODO |
[in,out] | sumtab | TODO |
Definition at line 953 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::layout, spmatrix_s::loc2glob, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::replicated, spmatrix_s::rowptr, spmatrix_s::values, and z_spm_oneinf_elt().
Referenced by z_spmOneInfNorm().
|
static |
Compute the one/inf norm of an spm CSR structure.
[in] | ntype | TODO |
[in] | spm | TODO |
[in,out] | sumtab | TODO |
Definition at line 1019 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::layout, spmatrix_s::loc2glob, spmatrix_s::mtxtype, spmatrix_s::n, spmatrix_s::replicated, spmatrix_s::rowptr, spmatrix_s::values, and z_spm_oneinf_elt().
Referenced by z_spmOneInfNorm().
|
static |
Compute the one/inf norm of an spm IJV structure.
[in] | ntype | TODO |
[in] | spm | TODO |
[in,out] | sumtab | TODO |
Definition at line 1085 of file z_spm_norm.c.
References spmatrix_s::baseval, spmatrix_s::colptr, spmatrix_s::dof, spmatrix_s::dofs, spmatrix_s::layout, spmatrix_s::mtxtype, spmatrix_s::nnz, spmatrix_s::rowptr, spmatrix_s::values, and z_spm_oneinf_elt().
Referenced by z_spmOneInfNorm().
|
static |
Compute the one/inf norm of the given spm structure given by the maximum row sum.
[in] | ntype | The type of norm to compute. |
[in] | spm | The spm from which the norm need to be computed. |
Definition at line 1149 of file z_spm_norm.c.
References spmatrix_s::clustnbr, spmatrix_s::comm, spmatrix_s::fmttype, spmatrix_s::gNexp, spmatrix_s::replicated, SpmCSC, SpmCSR, SpmIJV, z_spmOneInfNorm_csc(), z_spmOneInfNorm_csr(), and z_spmOneInfNorm_ijv().
Referenced by z_spmNorm().