Mentions légales du service

Skip to content
Snippets Groups Projects
Commit 0482c730 authored by Tony Delarue's avatar Tony Delarue Committed by Mathieu Faverge
Browse files

spmNorm is now distributed, and and a test function has been added

parent 55b02598
No related branches found
No related tags found
1 merge request!254 - Dist/SpmNorm
......@@ -11,7 +11,8 @@
* @author Xavier Lacoste
* @author Pierre Ramet
* @author Mathieu Faverge
* @date 2013-06-24
* @author Tony Delarue
* @date 2020-05-19
*
* @addtogroup spm
* @{
......@@ -505,33 +506,27 @@ double
spmNorm( spm_normtype_t ntype,
const spmatrix_t *spm )
{
spmatrix_t *spmtmp = (spmatrix_t*)spm;
double norm = -1.;
if ( spm->flttype == SpmPattern ) {
return SPM_ERR_BADPARAMETER;
return norm;
}
if ( spm->dof != 1 ) {
fprintf(stderr, "WARNING: spm expanded due to non implemented norm for non-expanded spm\n");
spmtmp = malloc( sizeof(spmatrix_t) );
spmExpand( spm, spmtmp );
}
switch (spm->flttype) {
case SpmFloat:
norm = (double)s_spmNorm( ntype, spmtmp );
norm = (double)s_spmNorm( ntype, spm );
break;
case SpmDouble:
norm = d_spmNorm( ntype, spmtmp );
norm = d_spmNorm( ntype, spm );
break;
case SpmComplex32:
norm = (double)c_spmNorm( ntype, spmtmp );
norm = (double)c_spmNorm( ntype, spm );
break;
case SpmComplex64:
norm = z_spmNorm( ntype, spmtmp );
norm = z_spmNorm( ntype, spm );
break;
case SpmPattern:
......@@ -539,10 +534,6 @@ spmNorm( spm_normtype_t ntype,
;
}
if ( spmtmp != spm ) {
spmExit( spmtmp );
free( spmtmp );
}
return norm;
}
......
This diff is collapsed.
......@@ -5,7 +5,7 @@
#
# @version 1.0.0
# @author Mathieu Faverge
# @date 2013-06-24
# @date 2020-05-19
#
###
include(RulesPrecisions)
......@@ -51,6 +51,7 @@ set (TESTS
if ( SPM_WITH_MPI )
list( APPEND TESTS
spm_scatter_gather_tests.c
spm_dist_norm_tests.c
)
endif()
......@@ -66,7 +67,8 @@ set( SPM_TESTS
set( SPM_DOF_TESTS
spm_dof_expand_tests spm_dof_norm_tests spm_dof_matvec_tests)
set( SPM_MPI_TESTS
spm_scatter_gather_tests )
spm_scatter_gather_tests
spm_dist_norm_tests )
# List of run types
set( RUNTYPE shm )
......
/**
*
* @file spm_dist_norm_tests.c
*
* Tests and validate the spm_norm routines.
*
* @copyright 2015-2017 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
* @version 6.0.0
* @author Mathieu Faverge
* @author Theophile Terraz
* @author Tony Delarue
* @date 2020-05-19
*
**/
#include <stdint.h>
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <string.h>
#include <assert.h>
#include <time.h>
#include <spm_tests.h>
#if !defined(SPM_WITH_MPI)
#error "This test should not be compiled in non distributed version"
#endif
int main (int argc, char **argv)
{
char *filename;
spmatrix_t original, *spm, *spmdist;
spm_driver_t driver;
int clustnbr = 1;
int clustnum = 0;
spm_mtxtype_t mtxtype;
spm_fmttype_t fmttype;
int baseval, distbycol, root;
int rc = SPM_SUCCESS;
int err = 0;
int dof, dofmax = 4;
int to_free = 0;
MPI_Init( &argc, &argv );
/**
* Get options from command line
*/
spmGetOptions( argc, argv,
&driver, &filename );
rc = spmReadDriver( driver, filename, &original );
free(filename);
if ( rc != SPM_SUCCESS ) {
fprintf(stderr, "ERROR: Could not read the file, stop the test !!!\n");
return EXIT_FAILURE;
}
MPI_Comm_size( MPI_COMM_WORLD, &clustnbr );
MPI_Comm_rank( MPI_COMM_WORLD, &clustnum );
if ( original.flttype == SpmPattern ) {
spmGenFakeValues( &original );
}
spmPrintInfo( &original, stdout );
printf(" -- SPM Norms Test --\n");
for( fmttype=SpmCSC; fmttype<=SpmIJV; fmttype++ ) {
distbycol = (fmttype == SpmCSR) ? 0 : 1;
if ( spmConvert( fmttype, &original ) != SPM_SUCCESS ) {
fprintf( stderr, "Issue to convert to %s format\n", fmtnames[fmttype] );
continue;
}
for( dof=-1; dof<2; dof++ )
{
if ( dof >= 0 ) {
spm = spmDofExtend( &original, dof, dofmax );
to_free = 1;
}
else {
spm = malloc(sizeof(spmatrix_t));
memcpy( spm, &original, sizeof(spmatrix_t) );
to_free = 0;
}
if ( spm == NULL ) {
fprintf( stderr, "Issue to extend the matrix\n" );
continue;
}
for( root=-1; root<clustnbr; root++ )
{
spmdist = spmScatter( spm, -1, NULL, distbycol, -1, MPI_COMM_WORLD );
if ( spmdist == NULL ) {
fprintf( stderr, "Failed to scatter the spm\n" );
err++;
continue;
}
for( baseval=0; baseval<2; baseval++ )
{
spmBase( spmdist, baseval );
for( mtxtype=SpmGeneral; mtxtype<=SpmHermitian; mtxtype++ )
{
if ( (mtxtype == SpmHermitian) &&
( ((original.flttype != SpmComplex64) &&
(original.flttype != SpmComplex32)) ||
(original.mtxtype != SpmHermitian) ) )
{
continue;
}
if ( (mtxtype != SpmGeneral) &&
(original.mtxtype == SpmGeneral) )
{
continue;
}
if ( spm ) {
spm->mtxtype = mtxtype;
}
spmdist->mtxtype = mtxtype;
if(clustnum == 0) {
printf( " Case: %s / %s / %d / %s / %d\n",
fltnames[spmdist->flttype],
fmtnames[spmdist->fmttype], baseval,
mtxnames[mtxtype - SpmGeneral], dof );
}
switch( spmdist->flttype ){
case SpmComplex64:
rc = z_spm_dist_norm_check( spm, spmdist );
break;
case SpmComplex32:
rc = c_spm_dist_norm_check( spm, spmdist );
break;
case SpmFloat:
rc = s_spm_dist_norm_check( spm, spmdist );
break;
case SpmDouble:
default:
rc = d_spm_dist_norm_check( spm, spmdist );
}
err = (rc != 0) ? err+1 : err;
}
}
spmExit( spmdist );
free( spmdist );
}
if ( spm != &original ) {
if( to_free ){
spmExit( spm );
}
free( spm );
}
}
}
MPI_Finalize();
if( err == 0 ) {
if(clustnum == 0) {
printf(" -- All tests PASSED --\n");
}
return EXIT_SUCCESS;
}
else
{
printf(" -- %d tests FAILED --\n", err);
return EXIT_FAILURE;
}
}
......@@ -9,7 +9,8 @@
*
* @version 1.0.0
* @author Mathieu Faverge
* @date 2013-06-24
* @author Tony Delarue
* @date 2020-05-19
*
**/
#ifndef _spm_tests_h_
......@@ -81,35 +82,42 @@ int core_sgeadd( spm_trans_t trans,
void z_spm_print_check( char *filename, const spmatrix_t *spm );
int z_spm_matvec_check( spm_trans_t trans, const spmatrix_t *spm );
int z_spm_norm_check( const spmatrix_t *spm );
int z_spm_dist_norm_check( const spmatrix_t *spm, const spmatrix_t *spmdist );
void c_spm_print_check( char *filename, const spmatrix_t *spm );
int c_spm_matvec_check( spm_trans_t trans, const spmatrix_t *spm );
int c_spm_norm_check( const spmatrix_t *spm );
int c_spm_dist_norm_check( const spmatrix_t *spm, const spmatrix_t *spmdist );
void d_spm_print_check( char *filename, const spmatrix_t *spm );
int d_spm_matvec_check( spm_trans_t trans, const spmatrix_t *spm );
int d_spm_norm_check( const spmatrix_t *spm );
int d_spm_dist_norm_check( const spmatrix_t *spm, const spmatrix_t *spmdist );
void s_spm_print_check( char *filename, const spmatrix_t *spm );
int s_spm_matvec_check( spm_trans_t trans, const spmatrix_t *spm );
int s_spm_norm_check( const spmatrix_t *spm );
int s_spm_dist_norm_check( const spmatrix_t *spm, const spmatrix_t *spmdist );
void p_spm_print_check( char *filename, const spmatrix_t *spm );
static inline int
spm_norm_print_result( double norms, double normd, double result )
spm_norm_print_result( double norms, double normd, double result, int clustnum )
{
int rc = 0;
if ( (result >= 0.) && (result < 1.) ) {
printf("SUCCESS !\n");
if(clustnum == 0) {
printf("SUCCESS !\n");
}
} else {
printf("FAILED !\n");
if(clustnum == 0) {
printf("FAILED !\n");
printf(" Nsparse = %e, Ndense = %e\n", norms, normd );
printf(" | Nsparse - Ndense | / Ndense = %e\n", result);
}
rc = 1;
}
printf(" Nsparse = %e, Ndense = %e\n", norms, normd );
printf(" | Nsparse - Ndense | / Ndense = %e\n", result);
return rc;
}
......
......@@ -10,7 +10,8 @@
* @version 6.0.0
* @author Mathieu Faverge
* @author Theophile Terraz
* @date 2015-01-01
* @author Tony Delarue
* @date 2020-05-19
*
* @precisions normal z -> c d s
*
......@@ -228,7 +229,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);
ret += spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result, 0 );
/**
* Test Norm Inf
......@@ -238,7 +239,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)) / nnz;
ret += spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result, 0 );
/**
* Test Norm One
......@@ -248,7 +249,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)) / nnz;
ret += spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result, 0 );
/**
* Test Norm Frobenius
......@@ -258,8 +259,79 @@ 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 / nnz;
ret += spm_norm_print_result( norms, normd, result );
ret += spm_norm_print_result( norms, normd, result, 0 );
free(A);
return ret;
}
/*------------------------------------------------------------------------
* Check the accuracy of the solution
*/
#if defined(SPM_WITH_MPI)
int
z_spm_dist_norm_check( const spmatrix_t *spm,
const spmatrix_t *spmdist )
{
double norms, normd;
double eps, result, nnz;
int ret = 0;
int clustnum = spmdist->clustnum;
eps = LAPACKE_dlamch_work('e');
nnz = spm->gnnzexp;
if ( spm->mtxtype != SpmGeneral ) {
nnz *= 2.;
}
/**
* Test Norm Max
*/
if(clustnum == 0) {
printf(" -- Test norm Max : ");
}
norms = spmNorm( SpmMaxNorm, spm );
normd = spmNorm( SpmMaxNorm, spmdist );
result = fabs(norms - normd) / (normd * eps);
ret += spm_norm_print_result( norms, normd, result, clustnum );
/**
* Test Norm Inf
*/
if(clustnum == 0) {
printf(" -- Test norm Inf : ");
}
norms = spmNorm( SpmInfNorm, spm );
normd = spmNorm( SpmInfNorm, spmdist );
result = fabs(norms - normd) / (normd * eps);
result = result * ((double)(spm->gNexp)) / nnz;
ret += spm_norm_print_result( norms, normd, result, clustnum );
/**
* Test Norm One
*/
if(clustnum == 0) {
printf(" -- Test norm One : ");
}
norms = spmNorm( SpmOneNorm, spm );
normd = spmNorm( SpmOneNorm, spmdist );
result = fabs(norms - normd) / (normd * eps);
result = result * ((double)(spm->gNexp)) / nnz;
ret += spm_norm_print_result( norms, normd, result, clustnum );
/**
* Test Norm Frobenius
*/
if(clustnum == 0) {
printf(" -- Test norm Frb : ");
}
norms = spmNorm( SpmFrobeniusNorm, spm );
normd = spmNorm( SpmFrobeniusNorm, spmdist );
result = fabs(norms - normd) / (normd * eps);
result = result / nnz;
ret += spm_norm_print_result( norms, normd, result, clustnum );
return ret;
}
#endif
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