Skip to content

Commit

Permalink
require an output location to be provided externally to eigsolveQuda …
Browse files Browse the repository at this point in the history
…in order to allow all requested eigenvalues to be returned, if desired
  • Loading branch information
kostrzewa committed Jul 19, 2023
1 parent 9c27ed5 commit a8b0a62
Show file tree
Hide file tree
Showing 3 changed files with 44 additions and 63 deletions.
54 changes: 28 additions & 26 deletions phmc.c
Original file line number Diff line number Diff line change
Expand Up @@ -211,7 +211,9 @@ void init_phmc() {
void phmc_compute_ev(const int trajectory_counter,
const int id,
matrix_mult_bi Qsq) {
double atime, etime, temp=0., temp2=0.;
double atime, etime;
double eval_min = 0.;
double eval_max = 0.;
int max_iter_ev, no_eigenvalues;
char buf[100];
char * phmcfilename = buf;
Expand All @@ -230,15 +232,15 @@ void phmc_compute_ev(const int trajectory_counter,
no_eigenvalues = 1;
if(mnl->external_eigsolver == QUDA_EIGSOLVER) {
#ifdef TM_USE_QUDA
temp = eigsolveQuda(no_eigenvalues, eigenvalue_precision, 1, 0, max_iter_ev, 0,
mnl->accprec, mnl->maxiter, mnl->eig_polydeg, mnl->eig_amin,
mnl->eig_amax, mnl->eig_n_kr, mnl->solver, g_relative_precision_flag,
1, // we only support even-odd here
mnl->solver_params.refinement_precision,
mnl->solver_params.sloppy_precision,
mnl->solver_params.compression_type, 0);
eigsolveQuda(&eval_min, no_eigenvalues, eigenvalue_precision, 1, 0, max_iter_ev, 0,
mnl->accprec, mnl->maxiter, mnl->eig_polydeg, mnl->eig_amin,
mnl->eig_amax, mnl->eig_n_kr, mnl->solver, g_relative_precision_flag,
1, // we only support even-odd here
mnl->solver_params.refinement_precision,
mnl->solver_params.sloppy_precision,
mnl->solver_params.compression_type, 0);
if( fabs(mnl->EVMax - 1) < 2*DBL_EPSILON ) {
temp = temp / mnl->StildeMax;
eval_min /= mnl->StildeMax;
}
#else
if(g_proc_id == 0) {
Expand All @@ -250,22 +252,22 @@ void phmc_compute_ev(const int trajectory_counter,
}
#endif
}else {
temp = eigenvalues_bi(&no_eigenvalues, max_iter_ev, eigenvalue_precision, 0, Qsq);
eval_min = eigenvalues_bi(&no_eigenvalues, max_iter_ev, eigenvalue_precision, 0, Qsq);
}


no_eigenvalues = 1;
if(mnl->external_eigsolver == QUDA_EIGSOLVER) {
#ifdef TM_USE_QUDA
temp2 = eigsolveQuda(no_eigenvalues, eigenvalue_precision, 1, 0, max_iter_ev, 1,
mnl->accprec, mnl->maxiter, mnl->eig_polydeg, mnl->eig_amin,
mnl->eig_amax, mnl->eig_n_kr, mnl->solver, g_relative_precision_flag,
1, // we only support even-odd here
mnl->solver_params.refinement_precision,
mnl->solver_params.sloppy_precision,
mnl->solver_params.compression_type, 0);
eigsolveQuda(&eval_max, no_eigenvalues, eigenvalue_precision, 1, 0, max_iter_ev, 1,
mnl->accprec, mnl->maxiter, mnl->eig_polydeg, mnl->eig_amin,
mnl->eig_amax, mnl->eig_n_kr, mnl->solver, g_relative_precision_flag,
1, // we only support even-odd here
mnl->solver_params.refinement_precision,
mnl->solver_params.sloppy_precision,
mnl->solver_params.compression_type, 0);
if( fabs(mnl->EVMax - 1.) < 2*DBL_EPSILON ) {
temp2 = temp2 / mnl->StildeMax;
eval_max /= mnl->StildeMax;
}
#else
if(g_proc_id == 0) {
Expand All @@ -277,26 +279,26 @@ void phmc_compute_ev(const int trajectory_counter,
}
#endif
}else {
temp2 = eigenvalues_bi(&no_eigenvalues, max_iter_ev, eigenvalue_precision, 1, Qsq);
eval_max = eigenvalues_bi(&no_eigenvalues, max_iter_ev, eigenvalue_precision, 1, Qsq);
}


if((g_proc_id == 0) && (g_debug_level > 1)) {
printf("# %s: lowest eigenvalue end of trajectory %d = %e\n",
mnl->name, trajectory_counter, temp);
mnl->name, trajectory_counter, eval_min);
printf("# %s: maximal eigenvalue end of trajectory %d = %e\n",
mnl->name, trajectory_counter, temp2);
mnl->name, trajectory_counter, eval_max);
}
if(g_proc_id == 0) {
if(temp2 > mnl->EVMax) {
fprintf(stderr, "\nWarning: largest eigenvalue for monomial %s larger than upper bound!\n\n", mnl->name);
if(eval_max > mnl->EVMax) {
fprintf(stderr, "\nWarning: largest eigenvalue for monomial %s: %.6f is larger than upper bound: %.6f\n\n", mnl->name, eval_max, mnl->EVMax);
}
if(temp < mnl->EVMin) {
fprintf(stderr, "\nWarning: smallest eigenvalue for monomial %s smaller than lower bound!\n\n", mnl->name);
if(eval_min < mnl->EVMin) {
fprintf(stderr, "\nWarning: smallest eigenvalue for monomial %s: %.6f is smaller than lower bound: %.6f\n\n", mnl->name, eval_min, mnl->EVMin);
}
countfile = fopen(phmcfilename, "a");
fprintf(countfile, "%.8d %1.5e %1.5e %1.5e %1.5e\n",
trajectory_counter, temp, temp2, mnl->EVMin, mnl->EVMax);
trajectory_counter, eval_min, eval_max, mnl->EVMin, mnl->EVMax);
fclose(countfile);
}
etime = gettime();
Expand Down
43 changes: 11 additions & 32 deletions quda_interface.c
Original file line number Diff line number Diff line change
Expand Up @@ -2611,11 +2611,11 @@ Interface function for Eigensolver on Quda
*********************************************************/


double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterations, int maxmin,
const double precision, const int max_iter, const int polydeg, const double amin,
const double amax, const int n_kr, const int solver_flag, const int rel_prec,
const int even_odd_flag, const SloppyPrecision refinement_precision,
SloppyPrecision sloppy_precision, CompressionType compression, const int oneFlavourFlag) {
void eigsolveQuda(double * evals, int n_evals, double tol, int blksize, int blkwise, int max_iterations, int maxmin,
const double precision, const int max_iter, const int polydeg, const double amin,
const double amax, const int n_kr, const int solver_flag, const int rel_prec,
const int even_odd_flag, const SloppyPrecision refinement_precision,
SloppyPrecision sloppy_precision, CompressionType compression, const int oneFlavourFlag) {

tm_stopwatch_push(&g_timers, __func__, "");

Expand Down Expand Up @@ -2648,17 +2648,6 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati
1 /*always QpQm*/);
}

// QUDA applies the MMdag operator, we need QpQm^{-1) in the end
// so we want QUDA to use the MdagM operator
inv_param.dagger = QUDA_DAG_YES;


_Complex double * eigenvls;
double returnvalue;

// allocate memory for eigenvalues
eigenvls = (_Complex double *)malloc((n)*sizeof(_Complex double));

// create new eig_param
eig_param = newQudaEigParam();

Expand All @@ -2675,6 +2664,8 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati
eig_param.invert_param->cuda_prec_eigensolver = inv_param.cuda_prec;
eig_param.invert_param->clover_cuda_prec_eigensolver = inv_param.clover_cuda_prec;

// for consistency with tmLQCD's own eigensolver we require a precision of at least
// 1e-14
if(tol < 1.e-14) {
eig_param.tol = 1.e-14;
eig_param.qr_tol = 1.e-14;
Expand All @@ -2683,8 +2674,6 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati
eig_param.qr_tol = tol;
}



if(blkwise == 1) {
eig_param.eig_type = QUDA_EIG_BLK_TR_LANCZOS;
eig_param.block_size = blksize;
Expand All @@ -2707,10 +2696,6 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati
eig_param.use_norm_op = QUDA_BOOLEAN_FALSE;
}

// BK: these defaults seem to work on a 32c64 ensemble
// at a relatively coarse lattice spacing for the eigenvalues
// of the twisted-clover ND operator with values of musigma / mudelta
// reproducing physical sea strange and charm quark masses
eig_param.use_poly_acc = (maxmin == 1) || (polydeg == 0) ? QUDA_BOOLEAN_FALSE : QUDA_BOOLEAN_TRUE;
eig_param.poly_deg = polydeg;
eig_param.a_min = amin;
Expand Down Expand Up @@ -2745,9 +2730,9 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati

/* The size of eigenvector search space and
* the number of required converged eigenvectors
* is both set to n */
eig_param.n_conv = n;
eig_param.n_ev = n;
* is both set to n_evals */
eig_param.n_conv = n_evals;
eig_param.n_ev = n_evals;
/* The size of the Krylov space is set to 96.
* From my understanding, QUDA automatically scales
* this search space, however more testing on this
Expand All @@ -2756,13 +2741,7 @@ double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterati

eig_param.max_restarts = max_iterations;

eigensolveQuda(NULL, eigenvls, &eig_param);

returnvalue = eigenvls[0];
free(eigenvls);
eigensolveQuda(NULL, evals, &eig_param);

tm_stopwatch_pop(&g_timers, 0, 1, "TM_QUDA");

return(returnvalue);

}
10 changes: 5 additions & 5 deletions quda_interface.h
Original file line number Diff line number Diff line change
Expand Up @@ -175,10 +175,10 @@ void compute_gauge_derivative_quda(monomial * const mnl, hamiltonian_field_t * c
void compute_WFlow_quda(const double eps ,const double tmax, const int traj, FILE* outfile);


double eigsolveQuda(int n, double tol, int blksize, int blkwise, int max_iterations, int maxmin,
const double precision, const int max_iter, const int polydeg, const double amin,
const double amax, const int n_kr, const int solver_flag, const int rel_prec,
const int even_odd_flag, const SloppyPrecision refinement_precision,
SloppyPrecision sloppy_precision, CompressionType compression, const int oneFlavourFlag);
void eigsolveQuda(double * evals, int n_evals, double tol, int blksize, int blkwise, int max_iterations, int maxmin,
const double precision, const int max_iter, const int polydeg, const double amin,
const double amax, const int n_kr, const int solver_flag, const int rel_prec,
const int even_odd_flag, const SloppyPrecision refinement_precision,
SloppyPrecision sloppy_precision, CompressionType compression, const int oneFlavourFlag);

#endif /* QUDA_INTERFACE_H_ */

0 comments on commit a8b0a62

Please sign in to comment.