Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

prepare for interface-breaking QUDA PR of feature/cheby-mg-setup #548

Open
wants to merge 7 commits into
base: master
Choose a base branch
from
14 changes: 12 additions & 2 deletions quda_dummy_types.h
Original file line number Diff line number Diff line change
Expand Up @@ -42,11 +42,11 @@ typedef enum QudaInverterType_s {
QUDA_CA_GCR_INVERTER = CA_GCR
} QudaInverterType;

typedef enum QudaCABasis_s {
typedef enum QudaPolynomialBasis_s {
QUDA_POWER_BASIS,
QUDA_CHEBYSHEV_BASIS,
QUDA_INVALID_BASIS = QUDA_INVALID_ENUM
} QudaCABasis;
} QudaPolynomialBasis;

typedef enum QudaEigSpectrumType_s {
QUDA_SPECTRUM_SR_EIG,
Expand Down Expand Up @@ -96,4 +96,14 @@ typedef struct QudaInverParam_s {
QudaPrecision clover_cuda_prec_eigensolver;
} QudaInvertParam;

typedef enum QudaNullVectorSetupType_s {
QUDA_SETUP_NULL_VECTOR_INVERSE_ITERATIONS,
QUDA_SETUP_NULL_VECTOR_CHEBYSHEV_FILTER,
QUDA_SETUP_NULL_VECTOR_EIGENVECTORS,
QUDA_SETUP_NULL_VECTOR_TEST_VECTORS,
QUDA_SETUP_NULL_VECTOR_RESTRICT_FINE,
QUDA_SETUP_NULL_VECTOR_FREE_FIELD,
QUDA_SETUP_NULL_VECTOR_INVALID = QUDA_INVALID_ENUM
} QudaNullVectorSetupType;

#endif
20 changes: 14 additions & 6 deletions quda_interface.c
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@
* 2018 Bartosz Kostrzewa, Ferenc Pittler
* 2019, 2020 Bartosz Kostrzewa
* 2021 Bartosz Kostrzewa, Marco Garofalo, Ferenc Pittler, Simone Bacchio
* 2022 Bartosz Kostrzewa, Marco Garofalo, Simone Romiti
*
* This file is part of tmLQCD.
*
Expand Down Expand Up @@ -1944,8 +1945,21 @@ void _setQudaMultigridParam(QudaMultigridParam* mg_param) {

} // for( dim=0 to dim=3 ) (space-time dimensions)

// set file i/o parameters to a default which imples that no null-vector I/O is performed
strcpy(mg_param->vec_infile[level], "");
strcpy(mg_param->vec_outfile[level], "");
mg_param->vec_store[level] = QUDA_BOOLEAN_NO;
mg_param->vec_load[level] = QUDA_BOOLEAN_NO;

mg_param->verbosity[level] = quda_input.mg_verbosity[level];
mg_param->precision_null[level] = QUDA_HALF_PRECISION;
// different setup types exist but inverse iterations are the default
// in the future we might explore also the Chebyshev filter with refinement to
// generate null vectors
mg_param->setup_type[level] = QUDA_SETUP_NULL_VECTOR_INVERSE_ITERATIONS;
mg_param->setup_restrict_remaining_type[level] = QUDA_SETUP_NULL_VECTOR_INVERSE_ITERATIONS;
mg_param->setup_maxiter_inverse_iterations_refinement[level] = 0;

mg_param->setup_inv_type[level] = quda_input.mg_setup_inv_type;
// Kate says: experimental, leave at 1 (will be used for bootstrap-style setup later)
mg_param->num_setup_iter[level] = 1;
Expand Down Expand Up @@ -2068,16 +2082,10 @@ void _setQudaMultigridParam(QudaMultigridParam* mg_param) {
// only coarsen the spin on the first restriction
mg_param->spin_block_size[0] = 2;

mg_param->compute_null_vector = QUDA_COMPUTE_NULL_VECTOR_YES;
mg_param->generate_all_levels = QUDA_BOOLEAN_YES;

mg_param->run_low_mode_check = quda_input.mg_run_low_mode_check;
mg_param->run_oblique_proj_check = quda_input.mg_run_oblique_proj_check;
mg_param->run_verify = quda_input.mg_run_verify;

// set file i/o parameters
strcpy(mg_param->vec_infile, "");
strcpy(mg_param->vec_outfile, "");
}

int invert_eo_degenerate_quda(spinor * const out,
Expand Down
134 changes: 69 additions & 65 deletions quda_types.h
Original file line number Diff line number Diff line change
Expand Up @@ -42,77 +42,81 @@ typedef enum tm_quda_ferm_bc_t {
} tm_quda_ferm_bc_t;

/* tm_QudaParams_t provides an interface between the tmLQCD input file and the
* available QUDA parameters. At the moment, only the fermionic bounday conditions
* and the MG parameters are exposed like this, but a further refactoring might
* turn this into a complete representation of the possible input parameters */
* available QUDA parameters for those where we don't
* a) rely on the defaults supplied by the various newQuda****Param() functions
* b) hard-code in quda_interface.c
* c) generate dynamically in quda_interface.c
* */
typedef struct tm_QudaParams_t {
int enable_device_memory_pool;
int enable_pinned_memory_pool;

tm_quda_ferm_bc_t fermionbc;

int pipeline;
double reliable_delta;
int gcrNkrylov;

int mg_n_level;
QudaVerbosity mg_verbosity[QUDA_MAX_MG_LEVEL];
int mg_n_vec[QUDA_MAX_MG_LEVEL];
int mg_blocksize[QUDA_MAX_MG_LEVEL][4];
double mg_mu_factor[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_setup_inv_type;
double mg_setup_2kappamu;
double mg_setup_tol[QUDA_MAX_MG_LEVEL];
int mg_setup_maxiter[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_coarse_solver_type[QUDA_MAX_MG_LEVEL];
int mg_coarse_solver_maxiter[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_tol[QUDA_MAX_MG_LEVEL];
int mg_nu_pre[QUDA_MAX_MG_LEVEL];
int mg_nu_post[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_smoother_type[QUDA_MAX_MG_LEVEL];
double mg_smoother_tol[QUDA_MAX_MG_LEVEL];
double mg_omega[QUDA_MAX_MG_LEVEL];
int mg_run_verify;
int mg_enable_size_three_blocks;
double mg_reuse_setup_mu_threshold;
double mg_reset_setup_mdu_threshold;
double mg_refresh_setup_mdu_threshold;

int mg_setup_maxiter_refresh[QUDA_MAX_MG_LEVEL];

int enable_device_memory_pool;
int enable_pinned_memory_pool;

tm_quda_ferm_bc_t fermionbc;

int pipeline;
double reliable_delta;
int gcrNkrylov;

int mg_n_level;
QudaVerbosity mg_verbosity[QUDA_MAX_MG_LEVEL];
int mg_n_vec[QUDA_MAX_MG_LEVEL];
int mg_blocksize[QUDA_MAX_MG_LEVEL][4];
double mg_mu_factor[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_setup_inv_type;
double mg_setup_2kappamu;
double mg_setup_tol[QUDA_MAX_MG_LEVEL];
int mg_setup_maxiter[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_coarse_solver_type[QUDA_MAX_MG_LEVEL];
int mg_coarse_solver_maxiter[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_tol[QUDA_MAX_MG_LEVEL];
int mg_nu_pre[QUDA_MAX_MG_LEVEL];
int mg_nu_post[QUDA_MAX_MG_LEVEL];
QudaInverterType mg_smoother_type[QUDA_MAX_MG_LEVEL];
double mg_smoother_tol[QUDA_MAX_MG_LEVEL];
double mg_omega[QUDA_MAX_MG_LEVEL];
int mg_run_verify;
int mg_enable_size_three_blocks;
double mg_reuse_setup_mu_threshold;
double mg_reset_setup_mdu_threshold;
double mg_refresh_setup_mdu_threshold;

int mg_setup_maxiter_refresh[QUDA_MAX_MG_LEVEL];

QudaNullVectorSetupType mg_setup_type[QUDA_MAX_MG_LEVEL];

// parameters related to communication-avoiding
// solvers
QudaCABasis mg_setup_ca_basis[QUDA_MAX_MG_LEVEL];
int mg_setup_ca_basis_size[QUDA_MAX_MG_LEVEL];
double mg_setup_ca_lambda_min[QUDA_MAX_MG_LEVEL];
double mg_setup_ca_lambda_max[QUDA_MAX_MG_LEVEL];
QudaPolynomialBasis mg_setup_ca_basis[QUDA_MAX_MG_LEVEL];
int mg_setup_ca_basis_size[QUDA_MAX_MG_LEVEL];
double mg_setup_ca_lambda_min[QUDA_MAX_MG_LEVEL];
double mg_setup_ca_lambda_max[QUDA_MAX_MG_LEVEL];

QudaCABasis mg_coarse_solver_ca_basis[QUDA_MAX_MG_LEVEL];
int mg_coarse_solver_ca_basis_size[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_ca_lambda_min[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_ca_lambda_max[QUDA_MAX_MG_LEVEL];
QudaPolynomialBasis mg_coarse_solver_ca_basis[QUDA_MAX_MG_LEVEL];
int mg_coarse_solver_ca_basis_size[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_ca_lambda_min[QUDA_MAX_MG_LEVEL];
double mg_coarse_solver_ca_lambda_max[QUDA_MAX_MG_LEVEL];

// parameters related to coarse grid deflation in the MG
int mg_use_eig_solver[QUDA_MAX_MG_LEVEL];
int mg_eig_preserve_deflation;
int mg_eig_nEv[QUDA_MAX_MG_LEVEL];
int mg_eig_nKr[QUDA_MAX_MG_LEVEL];
int mg_eig_require_convergence[QUDA_MAX_MG_LEVEL];
int mg_eig_check_interval[QUDA_MAX_MG_LEVEL];
int mg_eig_max_restarts[QUDA_MAX_MG_LEVEL];
double mg_eig_tol[QUDA_MAX_MG_LEVEL];
int mg_eig_use_poly_acc[QUDA_MAX_MG_LEVEL];
int mg_eig_poly_deg[QUDA_MAX_MG_LEVEL];
double mg_eig_amin[QUDA_MAX_MG_LEVEL];
double mg_eig_amax[QUDA_MAX_MG_LEVEL];
int mg_eig_use_normop[QUDA_MAX_MG_LEVEL];
int mg_eig_use_dagger[QUDA_MAX_MG_LEVEL];
QudaEigSpectrumType mg_eig_spectrum[QUDA_MAX_MG_LEVEL];
QudaEigType mg_eig_type[QUDA_MAX_MG_LEVEL];
int mg_coarse_guess;

int mg_run_low_mode_check;
int mg_run_oblique_proj_check;
int mg_use_eig_solver[QUDA_MAX_MG_LEVEL];
int mg_eig_preserve_deflation;
int mg_eig_nEv[QUDA_MAX_MG_LEVEL];
int mg_eig_nKr[QUDA_MAX_MG_LEVEL];
int mg_eig_require_convergence[QUDA_MAX_MG_LEVEL];
int mg_eig_check_interval[QUDA_MAX_MG_LEVEL];
int mg_eig_max_restarts[QUDA_MAX_MG_LEVEL];
double mg_eig_tol[QUDA_MAX_MG_LEVEL];
int mg_eig_use_poly_acc[QUDA_MAX_MG_LEVEL];
int mg_eig_poly_deg[QUDA_MAX_MG_LEVEL];
double mg_eig_amin[QUDA_MAX_MG_LEVEL];
double mg_eig_amax[QUDA_MAX_MG_LEVEL];
int mg_eig_use_normop[QUDA_MAX_MG_LEVEL];
int mg_eig_use_dagger[QUDA_MAX_MG_LEVEL];
QudaEigSpectrumType mg_eig_spectrum[QUDA_MAX_MG_LEVEL];
QudaEigType mg_eig_type[QUDA_MAX_MG_LEVEL];
int mg_coarse_guess;

int mg_run_low_mode_check;
int mg_run_oblique_proj_check;

} tm_QudaParams_t;

Expand Down
6 changes: 3 additions & 3 deletions read_input.l
Original file line number Diff line number Diff line change
Expand Up @@ -437,7 +437,7 @@ static inline double fltlist_next_token(int * const list_end){
free(input_copy);
}

static inline void parse_quda_mg_cabasis_par_array(char const * const input, QudaCABasis * par_array){
static inline void parse_quda_mg_polybasis_par_array(char const * const input, QudaPolynomialBasis * par_array){
char paramname[100];
char token[100];
char error_message[ERR_MSG_LEN];
Expand Down Expand Up @@ -1947,7 +1947,7 @@ static inline double fltlist_next_token(int * const list_end){
parse_quda_mg_dbl_par_array(yytext, &(quda_input.mg_eig_amax[0]));
}
{SPC}*MGSetupCABasisType{EQL}{STRLIST} {
parse_quda_mg_cabasis_par_array(yytext, &(quda_input.mg_setup_ca_basis[0]));
parse_quda_mg_polybasis_par_array(yytext, &(quda_input.mg_setup_ca_basis[0]));
}
{SPC}*MGSetupCABasisSize{EQL}{STRLIST} {
parse_quda_mg_int_par_array(yytext, &(quda_input.mg_setup_ca_basis_size[0]));
Expand All @@ -1959,7 +1959,7 @@ static inline double fltlist_next_token(int * const list_end){
parse_quda_mg_dbl_par_array(yytext, &(quda_input.mg_setup_ca_lambda_max[0]));
}
{SPC}*MGCoarseSolverCABasisType{EQL}{STRLIST} {
parse_quda_mg_cabasis_par_array(yytext, &(quda_input.mg_coarse_solver_ca_basis[0]));
parse_quda_mg_polybasis_par_array(yytext, &(quda_input.mg_coarse_solver_ca_basis[0]));
}
{SPC}*MGCoarseSolverCABasisSize{EQL}{STRLIST} {
parse_quda_mg_int_par_array(yytext, &(quda_input.mg_coarse_solver_ca_basis_size[0]));
Expand Down