Mentions légales du service

Skip to content
Snippets Groups Projects
Commit b5767493 authored by Mathieu Faverge's avatar Mathieu Faverge
Browse files

Fix issue in norm tests

parent 1fd19691
No related branches found
No related tags found
No related merge requests found
......@@ -57,7 +57,6 @@ sonar_spm:
expire_in: 1 week
paths:
- spm.lcov
- coverage/
- spm-coverage.xml
- spm-cppcheck.xml
- spm-rats.xml
......
......@@ -249,7 +249,7 @@ z_spmSymmetrize( spmatrix_t *spm )
spm_int_t lrow = oldcol[ j+1 ] - baseval;
int found = 0;
for (k = frow; (k < lrow); k++)
for (k = frow; k < lrow; k++)
{
if (i == (oldrow[k]-baseval))
{
......
......@@ -150,14 +150,14 @@ z_spmCSC2dense( const spmatrix_t *spm )
case SpmSymmetric:
for(j=0; j<spm->n; j++, colptr++)
{
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? (spm->dof * j) : dofs[j] - baseval;
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? spm->dof * j : dofs[j] - baseval;
for(k=colptr[0]; k<colptr[1]; k++, rowptr++)
{
i = (*rowptr - baseval);
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? (spm->dof * i) : dofs[i] - baseval;
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? spm->dof * i : dofs[i] - baseval;
for(jj=0; jj<dofj; jj++)
{
......@@ -174,14 +174,14 @@ z_spmCSC2dense( const spmatrix_t *spm )
default:
for(j=0; j<spm->n; j++, colptr++)
{
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? (spm->dof * j) : dofs[j] - baseval;
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? spm->dof * j : dofs[j] - baseval;
for(k=colptr[0]; k<colptr[1]; k++, rowptr++)
{
i = (*rowptr - baseval);
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? (spm->dof * i) : dofs[i] - baseval;
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? spm->dof * i : dofs[i] - baseval;
for(jj=0; jj<dofj; jj++)
{
......@@ -298,14 +298,14 @@ z_spmCSR2dense( const spmatrix_t *spm )
case SpmHermitian:
for(i=0; i<spm->n; i++, rowptr++)
{
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? (spm->dof * i) : dofs[i] - baseval;
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? spm->dof * i : dofs[i] - baseval;
for(k=rowptr[0]; k<rowptr[1]; k++, colptr++)
{
j = (*colptr - baseval);
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? (spm->dof * j) : dofs[j] - baseval;
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? spm->dof * j : dofs[j] - baseval;
for(jj=0; jj<dofj; jj++)
{
......@@ -328,14 +328,14 @@ z_spmCSR2dense( const spmatrix_t *spm )
case SpmSymmetric:
for(i=0; i<spm->n; i++, rowptr++)
{
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? (spm->dof * i) : dofs[i] - baseval;
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? spm->dof * i : dofs[i] - baseval;
for(k=rowptr[0]; k<rowptr[1]; k++, colptr++)
{
j = (*colptr - baseval);
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? (spm->dof * j) : dofs[j] - baseval;
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? spm->dof * j : dofs[j] - baseval;
for(jj=0; jj<dofj; jj++)
{
......@@ -352,14 +352,14 @@ z_spmCSR2dense( const spmatrix_t *spm )
default:
for(i=0; i<spm->n; i++, rowptr++)
{
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? (spm->dof * i) : dofs[i] - baseval;
dofi = ( spm->dof > 1 ) ? spm->dof : dofs[i+1] - dofs[i];
row = ( spm->dof > 1 ) ? spm->dof * i : dofs[i] - baseval;
for(k=rowptr[0]; k<rowptr[1]; k++, colptr++)
{
j = (*colptr - baseval);
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? (spm->dof * j) : dofs[j] - baseval;
dofj = ( spm->dof > 1 ) ? spm->dof : dofs[j+1] - dofs[j];
col = ( spm->dof > 1 ) ? spm->dof * j : dofs[j] - baseval;
for(jj=0; jj<dofj; jj++)
{
......
......@@ -118,7 +118,7 @@ static inline int
intsortcmp_iif( void ** const pbase, spm_int_t *p, spm_int_t *q ) {
spm_int_t *int1ptr = pbase[0];
spm_int_t *int2ptr = pbase[1];
return ( ( *p < *q ) || (( *p == *q ) && ( int2ptr[ p - int1ptr ] < int2ptr[ q - int1ptr ] )) );
return ( *p < *q ) || (( *p == *q ) && ( int2ptr[ p - int1ptr ] < int2ptr[ q - int1ptr ] ));
}
#define INTSORTCMP(p,q) intsortcmp_iif( pbase, (spm_int_t*)p, (spm_int_t*)q )
#include "integer_sort_mtypes.c"
......
......@@ -87,18 +87,21 @@ void s_spm_print_check( char *filename, const spmatrix_t *spm );
int s_spm_matvec_check( int trans, const spmatrix_t *spm );
int s_spm_norm_check( const spmatrix_t *spm );
static inline
void spm_norm_print_result( double norms, double normd, double result )
static inline int
spm_norm_print_result( double norms, double normd, double result )
{
int rc = 0;
if ( (result >= 0.) && (result < 1.) ) {
printf("SUCCESS !\n");
} else {
printf("FAILED !\n");
ret++;
rc=1;
}
printf(" Nsparse = %e, Ndense = %e\n", norms, normd );
printf(" | Nsparse - Ndense | / Ndense = %e\n", result);
return rc;
}
#endif /* _spm_tests_h_ */
......@@ -204,7 +204,7 @@ z_spm_norm_check( const spmatrix_t *spm )
norms = spmNorm( SpmMaxNorm, spm );
normd = LAPACKE_zlange( LAPACK_COL_MAJOR, 'M', spm->gNexp, spm->gNexp, A, spm->gNexp );
result = fabs(norms - normd) / (normd * eps);
spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result );
/**
* Test Norm Inf
......@@ -214,7 +214,7 @@ z_spm_norm_check( const spmatrix_t *spm )
normd = LAPACKE_zlange( LAPACK_COL_MAJOR, 'I', spm->gNexp, spm->gNexp, A, spm->gNexp );
result = fabs(norms - normd) / (normd * eps);
result = result * ((double)(spm->gNexp)) / ((double)(spm->gnnzexp));
spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result );
/**
* Test Norm One
......@@ -224,7 +224,7 @@ z_spm_norm_check( const spmatrix_t *spm )
normd = LAPACKE_zlange( LAPACK_COL_MAJOR, 'O', spm->gNexp, spm->gNexp, A, spm->gNexp );
result = fabs(norms - normd) / (normd * eps);
result = result * ((double)(spm->gNexp)) / ((double)(spm->gnnzexp));
spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result );
/**
* Test Norm Frobenius
......@@ -234,7 +234,7 @@ z_spm_norm_check( const spmatrix_t *spm )
normd = LAPACKE_zlange( LAPACK_COL_MAJOR, 'F', spm->gNexp, spm->gNexp, A, spm->gNexp );
result = fabs(norms - normd) / (normd * eps);
result = result / ((double)spm->gnnzexp);
spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result );
free(A);
return ret;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment