mld2p4-2:

config/pac.m4
 configure
 mlprec/impl/mld_dslud_interface.c
 mlprec/impl/mld_zslud_interface.c
 mlprec/mld_c_prec_type.f90
 mlprec/mld_d_prec_type.f90
 mlprec/mld_s_prec_type.f90
 mlprec/mld_z_prec_type.f90

Added coarse_solver tracker in prec_type.
Fixed configure and Xslud_interface for SuperLU_Dist version 5.
stopcriterion
Salvatore Filippone 8 years ago
parent cbe6426bef
commit 289fdc3b2a

@ -753,9 +753,26 @@ if test "x$pac_sludist_header_ok" == "xyes" ; then
LUstructInit(n, LUstruct);
}]])],
[ AC_MSG_RESULT([yes]); pac_sludist_version="4";],
[ AC_MSG_RESULT([no]); pac_sludist_version="2_3";])
[ AC_MSG_RESULT([no]); pac_sludist_version="3";])
AC_LANG_POP([C])
fi
if test "x$pac_sludist_version" == "x4" ; then
AC_MSG_CHECKING([for superlu_dist version 5])
AC_LANG_PUSH([C])
ac_cc=${MPICC-$CC}
AC_COMPILE_IFELSE(
[AC_LANG_SOURCE([[ #include "superlu_ddefs.h"
int testdslud()
{ superlu_dist_options_t options;
int n;
set_default_options_dist(&options);
}]])],
[ AC_MSG_RESULT([yes]); pac_sludist_version="5";],
[ AC_MSG_RESULT([no]); pac_sludist_version="4";])
AC_LANG_POP([C])
fi
fi
LIBS="$save_LIBS";
CPPFLAGS="$save_CPPFLAGS";
CC="$save_CC";

63
configure vendored

@ -12662,7 +12662,7 @@ else
sed 's/^/| /' conftest.$ac_ext >&5
{ $as_echo "$as_me:$LINENO: result: no" >&5
$as_echo "no" >&6; }; pac_sludist_version="2_3";
$as_echo "no" >&6; }; pac_sludist_version="3";
fi
rm -f core conftest.err conftest.$ac_objext conftest.$ac_ext
@ -12672,7 +12672,68 @@ ac_compile='$CC -c $CFLAGS $CPPFLAGS conftest.$ac_ext >&5'
ac_link='$CC -o conftest$ac_exeext $CFLAGS $CPPFLAGS $LDFLAGS conftest.$ac_ext $LIBS >&5'
ac_compiler_gnu=$ac_cv_c_compiler_gnu
if test "x$pac_sludist_version" == "x4" ; then
{ $as_echo "$as_me:$LINENO: checking for superlu_dist version 5" >&5
$as_echo_n "checking for superlu_dist version 5... " >&6; }
ac_ext=c
ac_cpp='$CPP $CPPFLAGS'
ac_compile='$CC -c $CFLAGS $CPPFLAGS conftest.$ac_ext >&5'
ac_link='$CC -o conftest$ac_exeext $CFLAGS $CPPFLAGS $LDFLAGS conftest.$ac_ext $LIBS >&5'
ac_compiler_gnu=$ac_cv_c_compiler_gnu
ac_cc=${MPICC-$CC}
cat >conftest.$ac_ext <<_ACEOF
/* confdefs.h. */
_ACEOF
cat confdefs.h >>conftest.$ac_ext
cat >>conftest.$ac_ext <<_ACEOF
/* end confdefs.h. */
#include "superlu_ddefs.h"
int testdslud()
{ superlu_dist_options_t options;
int n;
set_default_options_dist(&options);
}
_ACEOF
rm -f conftest.$ac_objext
if { (ac_try="$ac_compile"
case "(($ac_try" in
*\"* | *\`* | *\\*) ac_try_echo=\$ac_try;;
*) ac_try_echo=$ac_try;;
esac
eval ac_try_echo="\"\$as_me:$LINENO: $ac_try_echo\""
$as_echo "$ac_try_echo") >&5
(eval "$ac_compile") 2>conftest.er1
ac_status=$?
grep -v '^ *+' conftest.er1 >conftest.err
rm -f conftest.er1
cat conftest.err >&5
$as_echo "$as_me:$LINENO: \$? = $ac_status" >&5
(exit $ac_status); } && {
test -z "$ac_c_werror_flag" ||
test ! -s conftest.err
} && test -s conftest.$ac_objext; then
{ $as_echo "$as_me:$LINENO: result: yes" >&5
$as_echo "yes" >&6; }; pac_sludist_version="5";
else
$as_echo "$as_me: failed program was:" >&5
sed 's/^/| /' conftest.$ac_ext >&5
{ $as_echo "$as_me:$LINENO: result: no" >&5
$as_echo "no" >&6; }; pac_sludist_version="4";
fi
rm -f core conftest.err conftest.$ac_objext conftest.$ac_ext
ac_ext=c
ac_cpp='$CPP $CPPFLAGS'
ac_compile='$CC -c $CFLAGS $CPPFLAGS conftest.$ac_ext >&5'
ac_link='$CC -o conftest$ac_exeext $CFLAGS $CPPFLAGS $LDFLAGS conftest.$ac_ext $LIBS >&5'
ac_compiler_gnu=$ac_cv_c_compiler_gnu
fi
fi
LIBS="$save_LIBS";
CPPFLAGS="$save_CPPFLAGS";
CC="$save_CC";

@ -133,8 +133,13 @@ int mld_dsludist_fact(int n, int nl, int nnzl, int ffstr,
int i, panel_size, permc_spec, relax, info;
trans_t trans;
double drop_tol = 0.0, b[1], berr[1];
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;
int fst_row;
@ -157,7 +162,7 @@ int mld_dsludist_fact(int n, int nl, int nnzl, int ffstr,
ScalePermstruct = (ScalePermstruct_t *) SUPERLU_MALLOC(sizeof(ScalePermstruct_t));
LUstruct = (LUstruct_t *) SUPERLU_MALLOC(sizeof(LUstruct_t));
ScalePermstructInit(n,n, ScalePermstruct);
#if defined(SLUD_VERSION_4)
#if defined(SLUD_VERSION_4) || defined(SLUD_VERSION_5)
LUstructInit(n, LUstruct);
#elif defined(SLUD_VERSION_3)
LUstructInit(n,n, LUstruct);
@ -225,8 +230,13 @@ int mld_dsludist_solve(int itrans, int n, int nrhs,
trans_t trans;
double drop_tol = 0.0;
double *berr;
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;
@ -288,7 +298,7 @@ int mld_dsludist_free(void *f_factors)
* This routine can be called from Fortran.
*
* free all storage in the end
*
*
*/
#ifdef Have_SLUDist_
SuperMatrix *A;
@ -300,8 +310,13 @@ int mld_dsludist_free(void *f_factors)
trans_t trans;
double drop_tol = 0.0;
double *berr;
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;

@ -140,8 +140,13 @@ int mld_zsludist_fact(int n, int nl, int nnzl, int ffstr,
int i, panel_size, permc_spec, relax, info;
trans_t trans;
double drop_tol = 0.0,berr[1];
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;
int fst_row;
@ -164,7 +169,7 @@ int mld_zsludist_fact(int n, int nl, int nnzl, int ffstr,
ScalePermstruct = (ScalePermstruct_t *) SUPERLU_MALLOC(sizeof(ScalePermstruct_t));
LUstruct = (LUstruct_t *) SUPERLU_MALLOC(sizeof(LUstruct_t));
ScalePermstructInit(n,n, ScalePermstruct);
#if defined(SLUD_VERSION_4)
#if defined(SLUD_VERSION_4) || defined(SLUD_VERSION_5)
LUstructInit(n, LUstruct);
#elif defined(SLUD_VERSION_3)
LUstructInit(n,n, LUstruct);
@ -172,7 +177,6 @@ int mld_zsludist_fact(int n, int nl, int nnzl, int ffstr,
choke_on_me;
#endif
/* Set the default input options. */
set_default_options_dist(&options);
options.IterRefine=NO;
@ -238,8 +242,13 @@ int mld_zsludist_solve(int itrans, int n, int nrhs,
trans_t trans;
double drop_tol = 0.0;
double *berr;
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;
@ -313,8 +322,13 @@ int mld_zsludist_free(void *f_factors)
trans_t trans;
double drop_tol = 0.0;
double *berr;
mem_usage_t mem_usage;
#if defined(SLUD_VERSION_5)
superlu_dist_options_t options;
#elif defined(SLUD_VERSION_4)||defined(SLUD_VERSION_3)
superlu_options_t options;
#else
choke_on_me;
#endif
SuperLUStat_t stat;
factors_t *LUfactors;

@ -99,6 +99,16 @@ module mld_c_prec_type
! Number of outer sweeps. Sometimes 2 V-cycles may be better than 1 W-cycle.
!
integer(psb_ipk_) :: outer_sweeps = 1
!
! Coarse solver requires some tricky checks, and this needs we record the
! choice in the format given by the user, to keep track against what
! is put later in the multilevel array
!
integer(psb_ipk_) :: coarse_solver = -1
!
! The multilevel hierarchy
!
type(mld_c_onelev_type), allocatable :: precv(:)
contains
procedure, pass(prec) :: psb_c_apply2_vect => mld_c_apply2_vect

@ -99,6 +99,16 @@ module mld_d_prec_type
! Number of outer sweeps. Sometimes 2 V-cycles may be better than 1 W-cycle.
!
integer(psb_ipk_) :: outer_sweeps = 1
!
! Coarse solver requires some tricky checks, and this needs we record the
! choice in the format given by the user, to keep track against what
! is put later in the multilevel array
!
integer(psb_ipk_) :: coarse_solver = -1
!
! The multilevel hierarchy
!
type(mld_d_onelev_type), allocatable :: precv(:)
contains
procedure, pass(prec) :: psb_d_apply2_vect => mld_d_apply2_vect

@ -99,6 +99,16 @@ module mld_s_prec_type
! Number of outer sweeps. Sometimes 2 V-cycles may be better than 1 W-cycle.
!
integer(psb_ipk_) :: outer_sweeps = 1
!
! Coarse solver requires some tricky checks, and this needs we record the
! choice in the format given by the user, to keep track against what
! is put later in the multilevel array
!
integer(psb_ipk_) :: coarse_solver = -1
!
! The multilevel hierarchy
!
type(mld_s_onelev_type), allocatable :: precv(:)
contains
procedure, pass(prec) :: psb_s_apply2_vect => mld_s_apply2_vect

@ -99,6 +99,16 @@ module mld_z_prec_type
! Number of outer sweeps. Sometimes 2 V-cycles may be better than 1 W-cycle.
!
integer(psb_ipk_) :: outer_sweeps = 1
!
! Coarse solver requires some tricky checks, and this needs we record the
! choice in the format given by the user, to keep track against what
! is put later in the multilevel array
!
integer(psb_ipk_) :: coarse_solver = -1
!
! The multilevel hierarchy
!
type(mld_z_onelev_type), allocatable :: precv(:)
contains
procedure, pass(prec) :: psb_z_apply2_vect => mld_z_apply2_vect

Loading…
Cancel
Save