#include <QMCFlags.h>
Public Member Functions | |
QMCFlags () | |
Creates an instance of this class. | |
void | read_flags (string InFileName) |
Load this object's state from a QMC input file and initialize the object. | |
void | set_filenames (string runfile) |
From the input file name, generate other file names that will be used during the calculation. | |
bool | checkFlags () |
Check for some obvious flaws with parameter combinations. | |
int | getNumGPUWalkers () |
If we are going to allow the GPU to calculate a subset of walkers out of the walkers_per_pass batch, then we need to know how many it is going to process. | |
Public Attributes | |
vector< long > | programmersLongs |
This facilitates easy testing of parameters through the input file. | |
string | base_file_name |
Base file name for the calculation. | |
string | input_file_name |
Input file name for the calculation. | |
string | output_file_name |
Output file name for the calculation. | |
string | checkin_file_name |
Name of the checkpoint file to read in. | |
string | checkout_file_name |
Name of the checkpoint file to write out. | |
string | restart_file_name |
Restart input file name for the calculation. | |
string | results_file_name |
Output results file name for the calculation. | |
string | config_file_name |
Saved walker file name for the calculation. | |
string | density_file_name |
Basis function density file name for the calculation. | |
string | force_file_name |
Storage for the calculated forces. | |
string | walker_initialization_method |
Method for initializing the walkers for the calculation. | |
int | walker_initialization_combinations |
Number of walkers to generate when initializing the walkers. | |
string | run_type |
Type of QMC calculation. | |
int | one_e_per_iter |
Whether to update one electron at a time (1) or all at once (0). | |
int | use_surfer |
Using QMCSurfer instead of QMCWalker. | |
int | set_debug |
Minimize the amount of output. | |
int | use_jastrow |
These are debugging flags. | |
int | detailed_energies |
string | temp_dir |
Scratch directory for the calculation. | |
string | parallelization_method |
Parallelization algorithm to use. | |
long | iseed |
Seed used in generating pseudo-random numbers for the calculation. | |
string | nuclear_derivatives |
Implemented here is the method described in: S. | |
string | sampling_method |
Diffusion Green's function to be used for the calculation. | |
string | QF_modification_type |
Should the "quantum force" for the calculation be modified. | |
string | energy_modification_type |
Should the "local energy" for the calculation be modified. | |
double | umrigar93_equalelectrons_parameter |
If Umrigar's 1993 QMC algorithm (J. | |
string | walker_reweighting_method |
Algorithm used to reweight the walkers after a time step. | |
double | trail_eps2 |
Trail's eps^2 value for his weighting method found in Alternative sampling for variational quantum Monte Carlo Phys Rev E 77, 016704 (2008). | |
double | e_0 |
The e_0 value used in: Alternative sampling for variational quantum Monte Carlo Phys Rev E 77, 016704 (2008). | |
string | branching_method |
Method for branching walkers. | |
double | branching_threshold |
If a walker's weight exceeds this threshold, it will branch when the branching method allows the walkers to have fractional weights. | |
double | fusion_threshold |
If a walker's weight is less than this threshold, it will fuse with another low weight walker when the branching method allows the walkers to have fractional weights. | |
int | synchronize_dmc_ensemble |
1 if the global results are to be collected and broadcast to all nodes and used to calculate the estimated energy, trial energy, and effective time step. | |
int | synchronize_dmc_ensemble_interval |
The interval at which the DMC ensemble is to be synchronized. | |
int | old_walker_acceptance_parameter |
Parameter which can be used to prevent persistent configurations from interfering with a calculation. | |
int | warn_verbosity |
Warning verbosity level. | |
double | rel_cutoff |
Controls. | |
int | limit_branching |
Branching recommendations. | |
double | dt |
Time step being used for the calculation. | |
double | accel_delta |
These two parameters are used by Umrigar's metropolis acceleration scheme, labeled here umrigar93_accelerated_sampling. | |
double | accel_tm |
double | dt_run |
Time step to be used for the non-equilibration portion of the calculation. | |
double | dt_effective |
Last updated value of the effective time step. | |
double | desired_convergence |
Calculation terminates if the standard deviation of the energy falls below this threshold. | |
unsigned long | max_time_steps |
Maximum number of time steps the calculation will run for. | |
unsigned long | original_max_time_steps |
double | max_time |
Maximum amount of time the calculation will run for. | |
long | walkers_per_pass |
This parameter is used in QMCRun to process several walkers simultaneously. | |
vector< int > | future_walking |
Implemented here is the future walking method described in: J. | |
long | gpu_walkers_per_pass |
long | number_of_walkers |
Current number of walkers. | |
long | number_of_walkers_initial |
Initial number of walkers. | |
int | use_basis_function_interpolation |
Use interpolation of the radial basis functions in an attempt to speed up the calculation. | |
int | number_basis_function_interpolation_grid_points |
Number of grid points to use in interpolating the radial basis functions. | |
double | basis_function_interpolation_first_point |
Smallest point to use in interpolating the radial basis functions. | |
int | psuedo_gridLevel |
This is the index of the Lebedev-Laikov grid to use for the run on each nucleus. | |
double | psuedo_cutoff |
If the local part of the psuedopotential is sufficiently small, then we can forgo the calculation of the nonlocal part. | |
int | use_psuedopotential |
This is an internal flag, not an input flag. | |
int | output_interval |
Number of time steps taken between producing output. | |
int | zero_out_checkpoint_statistics |
1 if the accumulated statistics in the input checkpoint file are to be discarded, and 0 otherwise. | |
unsigned long | equilibration_steps |
Number of time steps to take for the calculation to equilibrate. | |
double | dt_equilibration |
Initial time step used during the equilibration of the calculation. | |
string | equilibration_function |
Procedure used to modify the time step while the calculation is equilibrating. | |
unsigned int | CKAnnealingEquilibration1_parameter |
Number of time steps to use the CKAnnealingEquilibration1 before just equilibrating with the time step of the non-equilibration portion of the calculation. | |
int | use_equilibration_array |
1 if QMCEquilibrationArray is to be used. | |
int | mpireduce_interval |
Number of time steps taken on the root processor before the results are collected from all other processors. | |
int | mpipoll_interval |
Number of time steps taken on the worker processors before checking for commands from the root processor. | |
int | checkpoint_interval |
Number of time steps taken before this processor writes a checkpoint file. | |
int | checkpoint |
1 if the calculation is to checkpoint, and 0 otherwise. | |
int | use_available_checkpoints |
Use existing checkpoints when performing this calculation. | |
int | checkpoint_energy_only |
To make them smaller, write only the energy data to the checkpoint files. | |
int | print_transient_properties_interval |
Number of time steps between printing the global properties as the calculation progresses. | |
int | print_transient_properties |
1 if the global properties of the calculation are output as the calculation progresses, and 0 otherwise. | |
int | Natoms |
Number of atoms in the system. | |
int | charge |
Charge of the system in atomic units. | |
int | Norbitals |
Number of orbitals in the trial wavefunction. | |
int | Nbasisfunc |
Number of basis functions used to create the orbitals for the trial wavefunction. | |
int | Ndeterminants |
Number of determinants in the SCF part of the wavefunction. | |
string | trial_function_type |
Restricted or unrestricted trial functions can be used. | |
int | calculate_bf_density |
1 if the basis function density is to be calculated, 0 if not. | |
double | energy_trial |
Trial energy used for branching QMC calculations. | |
double | energy_estimated |
Last estimate of the energy of the system. | |
double | energy_estimated_original |
Original estimate of the energy of the system read from the input file. | |
int | correct_population_size_bias |
1 if the correction for the population size bias is used, and 0 otherwise. | |
int | print_configs |
1 if walkers from the calculation are written out to file, and 0 otherwise. | |
int | print_config_frequency |
Number of time steps between writing walkers out to file. | |
int | optimize_Psi |
1 if the wavefunction is to be optimized during a VMC wavefunction, and 0 otherwise. | |
int | max_optimize_Psi_steps |
Number of vmc-optimize iterations used when optimizing the wavefunction. | |
double | optimize_Psi_barrier_parameter |
Objective function parameter used in forming a barrier around a valid region of parameter space. | |
int | calculate_Derivatives |
This is an internal flag. | |
int | optimize_EE_Jastrows |
These next flags control which part of the wavefunction we're optimizing. | |
int | optimize_EN_Jastrows |
int | optimize_NEE_Jastrows |
int | optimize_CI |
int | optimize_Orbitals |
int | optimize_L |
Some Jastrows (e.g. | |
int | use_three_body_jastrow |
int | reproduce_NE_with_NEE_jastrow |
If this flag is true, we allow the terms in the three body Jastrow that overlap the NE Jastrow to be nonzero. | |
int | reproduce_EE_with_NEE_jastrow |
If this flag is true, we allow the terms in the three body Jastrow that overlap the EE Jastrow to be nonzero. | |
double | a_diag |
This controls the initial value used by QMCLineSearch for modifying the hessian. | |
double | ksi |
A parameter used by the Linearize step length method, must be between 0 and 1. | |
string | optimize_Psi_criteria |
Objective function to optimize when optimizing a VMC wavefunction. | |
string | optimize_Psi_method |
Numerical optimization algorithm to use when optimizing a VMC wavefunction. | |
string | numerical_derivative_surface |
Objective function to use when calculating derivatives for optimizing a VMC wavefunction. | |
int | optimization_max_iterations |
Maximum number of iterations to use when optimizing a VMC wavefunction with one set of correlated-sampling configurations. | |
double | optimization_error_tolerance |
Tolerance used to determine when a numerical optimization with one set of correlated-sampling configurations has converged. | |
int | ck_genetic_algorithm_1_population_size |
Population size used when the CKGeneticAlgorithm1 is used to optimize a VMC wavefunction. | |
double | ck_genetic_algorithm_1_mutation_rate |
Mutation rate used when the CKGeneticAlgorithm1 is used to optimize a VMC wavefunction. | |
double | ck_genetic_algorithm_1_initial_distribution_deviation |
Standard deviation used when initializing the CKGeneticAlgorithm1. | |
string | line_search_step_length |
Algorithm used to determine the step length used when a line search algorithm is used to optimize a VMC wavefunction. | |
double | singularity_penalty_function_parameter |
Parameter used for the logarithmic barrier when D.R. | |
int | link_Jastrow_parameters |
1 if up and down electrons use the same Jastrow parameters for their electron-nucleus terms, and 0 otherwise. | |
int | link_NEE_Jastrows |
int | link_Orbital_parameters |
1 if up and down electrons use the same orbital parameters. | |
int | constrain_Orbital_zeros |
int | constrain_Orbital_same |
int | link_Determinant_parameters |
Just because the determinants shared orbitals when they were defined using the input file, that doesn't mean we can't optimize them independently. | |
int | replace_electron_nucleus_cusps |
1 if Gaussian orbitals are to be replaced by exponentials that satisfy the electron nucleus cusp consitions near nuclei and 0 otherwise. | |
string | energy_cutoff_type |
1 if we'll use the energy cutoff recommended by DePasquale, Rothstein, Vrbik JCP, 89, 3629, 1988 0 otherwise | |
int | equilibrate_every_opt_step |
1 if every VMC calculation used to generated correlated-sampling configurations for an optimization equilibrates for the specified equilibration time, and 0 otherwise. | |
int | equilibrate_first_opt_step |
1 if the first VMC calculation used to generated correlated-sampling configurations for an optimization equilibrates for the specified equilibration time, and 0 otherwise. | |
double | population_control_parameter |
Parameter used in controlling the number of walkers or total weight of all walkers during a branching QMC calculation. | |
int | write_all_energies_out |
1 if the local energy for every walker for every time step is to be written out, and 0 otherwise. | |
int | write_electron_densities |
1 if the distance between each pair of electrons for every time step is to be collected and written out in a histogram, and 0 otherwise. | |
double | max_pair_distance |
The maximum distance, in atomic units, of the histograms for the electron densities. | |
int | writePllxCorrelationDiagram |
1 if the x coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | pllxCorrelationDiagramMin |
The minimum value for the x coordinate of the 2D parallel spin correlation diagram. | |
double | pllxCorrelationDiagramMax |
The maximum value for the x coordinate of the 2D parallel spin correlation diagram. | |
int | writePllyCorrelationDiagram |
1 if the y coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | pllyCorrelationDiagramMin |
The minimum value for the y coordinate of the 2D parallel spin correlation diagram. | |
double | pllyCorrelationDiagramMax |
The maximum value for the y coordinate of the 2D parallel spin correlation diagram. | |
int | writePllzCorrelationDiagram |
1 if the z coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | pllzCorrelationDiagramMin |
The minimum value for the z coordinate of the 2D parallel spin correlation diagram. | |
double | pllzCorrelationDiagramMax |
The maximum value for the z coordinate of the 2D parallel spin correlation diagram. | |
int | writeOppxCorrelationDiagram |
1 if the x coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | oppxCorrelationDiagramMin |
The minimum value for the x coordinate of the 2D opposite spin correlation diagram. | |
double | oppxCorrelationDiagramMax |
The maximum value for the x coordinate of the 2D opposite spin correlation diagram. | |
int | writeOppyCorrelationDiagram |
1 if the y coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | oppyCorrelationDiagramMin |
The minimum value for the y coordinate of the 2D opposite spin correlation diagram. | |
double | oppyCorrelationDiagramMax |
The maximum value for the y coordinate of the 2D opposite spin correlation diagram. | |
int | writeOppzCorrelationDiagram |
1 if the z coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise. | |
double | oppzCorrelationDiagramMin |
The minimum value for the z coordinate of the 2D opposite spin correlation diagram. | |
double | oppzCorrelationDiagramMax |
The maximum value for the z coordinate of the 2D opposite spin correlation diagram. | |
string | chip_and_mike_are_cool |
Are Chip and Mike cool? Answer: Yea Baby! | |
int | my_rank |
MPI rank of this processor. | |
int | nprocs |
Number of processors used in this calculation. | |
int | use_hf_potential |
Hartree-Fock calculations using QMC, aimed at producing basis independent reference wavefunctions. | |
int | hf_num_average |
Number of electron samples to average v_eff over. | |
int | lock_trial_energy |
Option to keep trial energy fixed to its original value. | |
Friends | |
ostream & | operator<< (ostream &strm, QMCFlags &flags) |
Write this object's state out to a stream. |
Definition at line 30 of file QMCFlags.h.
QMCFlags::QMCFlags | ( | ) |
void QMCFlags::read_flags | ( | string | InFileName | ) |
Load this object's state from a QMC input file and initialize the object.
InFileName | QMC input file to load this object's state from. |
Definition at line 24 of file QMCFlags.cpp.
References a_diag, accel_delta, accel_tm, basis_function_interpolation_first_point, branching_method, branching_threshold, calculate_bf_density, calculate_Derivatives, charge, checkFlags(), checkin_file_name, checkout_file_name, checkpoint, checkpoint_energy_only, checkpoint_interval, chip_and_mike_are_cool, ck_genetic_algorithm_1_initial_distribution_deviation, ck_genetic_algorithm_1_mutation_rate, ck_genetic_algorithm_1_population_size, CKAnnealingEquilibration1_parameter, constrain_Orbital_same, constrain_Orbital_zeros, correct_population_size_bias, desired_convergence, detailed_energies, dt, dt_effective, dt_equilibration, dt_run, energy_cutoff_type, energy_estimated, energy_estimated_original, energy_modification_type, energy_trial, equilibrate_every_opt_step, equilibrate_first_opt_step, equilibration_function, equilibration_steps, fusion_threshold, future_walking, gpu_walkers_per_pass, hf_num_average, iseed, ksi, limit_branching, line_search_step_length, link_Determinant_parameters, link_Jastrow_parameters, link_NEE_Jastrows, link_Orbital_parameters, lock_trial_energy, max_optimize_Psi_steps, max_pair_distance, max_time, max_time_steps, mpipoll_interval, mpireduce_interval, Natoms, Nbasisfunc, Ndeterminants, Norbitals, nuclear_derivatives, number_basis_function_interpolation_grid_points, number_of_walkers, number_of_walkers_initial, numerical_derivative_surface, old_walker_acceptance_parameter, one_e_per_iter, oppxCorrelationDiagramMax, oppxCorrelationDiagramMin, oppyCorrelationDiagramMax, oppyCorrelationDiagramMin, oppzCorrelationDiagramMax, oppzCorrelationDiagramMin, optimization_error_tolerance, optimization_max_iterations, optimize_CI, optimize_EE_Jastrows, optimize_EN_Jastrows, optimize_L, optimize_NEE_Jastrows, optimize_Orbitals, optimize_Psi, optimize_Psi_barrier_parameter, optimize_Psi_criteria, optimize_Psi_method, original_max_time_steps, output_interval, parallelization_method, pllxCorrelationDiagramMax, pllxCorrelationDiagramMin, pllyCorrelationDiagramMax, pllyCorrelationDiagramMin, pllzCorrelationDiagramMax, pllzCorrelationDiagramMin, population_control_parameter, print_config_frequency, print_configs, print_transient_properties, print_transient_properties_interval, programmersLongs, psuedo_cutoff, psuedo_gridLevel, QF_modification_type, rel_cutoff, replace_electron_nucleus_cusps, reproduce_EE_with_NEE_jastrow, reproduce_NE_with_NEE_jastrow, run_type, sampling_method, set_debug, set_filenames(), singularity_penalty_function_parameter, synchronize_dmc_ensemble, synchronize_dmc_ensemble_interval, temp_dir, trial_function_type, umrigar93_equalelectrons_parameter, use_available_checkpoints, use_basis_function_interpolation, use_equilibration_array, use_hf_potential, use_jastrow, use_psuedopotential, use_surfer, use_three_body_jastrow, walker_initialization_combinations, walker_initialization_method, walker_reweighting_method, walkers_per_pass, warn_verbosity, write_all_energies_out, write_electron_densities, writeOppxCorrelationDiagram, writeOppyCorrelationDiagram, writeOppzCorrelationDiagram, writePllxCorrelationDiagram, writePllyCorrelationDiagram, writePllzCorrelationDiagram, and zero_out_checkpoint_statistics.
Referenced by QMCInput::read().
void QMCFlags::set_filenames | ( | string | runfile | ) |
From the input file name, generate other file names that will be used during the calculation.
(e.g. base_file_name)
runfile | QMC input file name |
Definition at line 883 of file QMCFlags.cpp.
References base_file_name, checkin_file_name, checkout_file_name, config_file_name, density_file_name, force_file_name, input_file_name, my_rank, output_file_name, restart_file_name, results_file_name, and temp_dir.
Referenced by read_flags().
bool QMCFlags::checkFlags | ( | ) |
Check for some obvious flaws with parameter combinations.
If something is irreparably wrong, then it will return false.
Otherwise, it will try to correct any problems automatically and return true;
In allowing a harmonic oscillator model, we use some of the regular molecular parameters to help define our HO system.
Definition at line 1201 of file QMCFlags.cpp.
References accel_delta, accel_tm, calculate_Derivatives, charge, checkout_file_name, checkpoint, chip_and_mike_are_cool, dt, dt_equilibration, equilibration_steps, QMCLineSearchStepLengthSelectionFactory::factory(), QMCInput::flags, getNumGPUWalkers(), globalInput, gpu_walkers_per_pass, ksi, line_search_step_length, link_Orbital_parameters, max_time_steps, Natoms, Nbasisfunc, nuclear_derivatives, number_of_walkers, one_e_per_iter, optimization_max_iterations, optimize_Psi, optimize_Psi_criteria, optimize_Psi_method, print_configs, psuedo_gridLevel, rel_cutoff, replace_electron_nucleus_cusps, run_type, sampling_method, temp_dir, trial_function_type, walker_initialization_method, and walkers_per_pass.
Referenced by read_flags().
int QMCFlags::getNumGPUWalkers | ( | ) | [inline] |
If we are going to allow the GPU to calculate a subset of walkers out of the walkers_per_pass batch, then we need to know how many it is going to process.
If we never set this parameter in the input file, then set it to be the entire walkers_per_pass batch.
If QMC_GPU is not set (we aren't using the GPU, then return ZERO, since the GPU will be calculating ZERO of the results).
Definition at line 1111 of file QMCFlags.h.
References gpu_walkers_per_pass, and walkers_per_pass.
Referenced by QMCSlater::allocate(), checkFlags(), QMCSCFJastrow::evaluate(), QMCSlater::initialize(), QMCJastrow::initialize(), and QMCSlater::update_Ds().
ostream& operator<< | ( | ostream & | strm, | |
QMCFlags & | flags | |||
) | [friend] |
Write this object's state out to a stream.
The same format is used as in the QMC input file.
Definition at line 968 of file QMCFlags.cpp.
vector<long> QMCFlags::programmersLongs |
This facilitates easy testing of parameters through the input file.
Enter as many as desired, separated by a space.
Definition at line 42 of file QMCFlags.h.
Referenced by read_flags().
string QMCFlags::base_file_name |
Base file name for the calculation.
This is the file name without any extension.
Definition at line 48 of file QMCFlags.h.
Referenced by qmcbeaver(), QMCManager::run(), set_filenames(), QMCManager::writeElectronDensityHistograms(), and QMCManager::writeTransientProperties().
string QMCFlags::input_file_name |
Input file name for the calculation.
This is the base file name with the ".ckmf" extension.
Definition at line 54 of file QMCFlags.h.
Referenced by QMCPropertyArrays::QMCPropertyArrays(), set_filenames(), and QMCManager::writeRestart().
string QMCFlags::output_file_name |
Output file name for the calculation.
This is the base file name with the ".qmc" extension.
Definition at line 60 of file QMCFlags.h.
Referenced by QMCManager::initializeOutputs(), and set_filenames().
string QMCFlags::checkin_file_name |
Name of the checkpoint file to read in.
The checkpoint files will automatically have a suffix added to them according to parallel node rank. It will search both the current directory for the checkpoint file, as well as in a directory also called checkin_file_name.
There is some use in permitting a different input than output since several runs might want to use the same input, and you don't want to overwrite the original checkpoints. For example, the checkpoint represents a (probably) equilibrated distribution.
If this flag is not specified, it will assume the name of the ckmf file.
Definition at line 77 of file QMCFlags.h.
Referenced by QMCManager::initializeCalculationState(), read_flags(), and set_filenames().
string QMCFlags::checkout_file_name |
Name of the checkpoint file to write out.
It will try to put the checkpoint files into a directory with the name checkout_file_name.
If this flag is not specified, it will assume the name of the ckmf file.
Definition at line 87 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), qmcbeaver(), read_flags(), set_filenames(), and QMCManager::writeCheckpoint().
string QMCFlags::restart_file_name |
Restart input file name for the calculation.
This is the base file name with the ".01.ckmf" extension.
Definition at line 93 of file QMCFlags.h.
Referenced by set_filenames(), and QMCManager::writeRestart().
string QMCFlags::results_file_name |
Output results file name for the calculation.
This is the base file name with the ".rslts" extension.
Definition at line 99 of file QMCFlags.h.
Referenced by QMCManager::initializeOutputs(), and set_filenames().
string QMCFlags::config_file_name |
Saved walker file name for the calculation.
This is the base file name with the ".<processor>.cfgs" extension.
Definition at line 105 of file QMCFlags.h.
Referenced by QMCReadAndEvaluateConfigs::locally_CalculateProperties(), QMCInput::openConfigFile(), and set_filenames().
string QMCFlags::density_file_name |
Basis function density file name for the calculation.
This is the base file name with the ".density" extension.
Definition at line 111 of file QMCFlags.h.
Referenced by set_filenames(), and QMCManager::writeBFDensity().
string QMCFlags::force_file_name |
Storage for the calculated forces.
This is the base file name with the ".force" extension.
Definition at line 117 of file QMCFlags.h.
Referenced by set_filenames(), and QMCManager::writeForces().
Method for initializing the walkers for the calculation.
Definition at line 123 of file QMCFlags.h.
Referenced by checkFlags(), QMCWalker::initializeWalkerPosition(), operator<<(), QMCRun::randomlyInitializeWalkers(), and read_flags().
Number of walkers to generate when initializing the walkers.
The desired number of walkers will be chosen from the generated walkers.
Definition at line 129 of file QMCFlags.h.
Referenced by QMCMikesBetterWalkerInitialization::initializeBunchOfWalkersPosition(), QMCMikesBetterWalkerInitialization::initializeWalkerPosition(), operator<<(), and read_flags().
string QMCFlags::run_type |
Type of QMC calculation.
For exampe: variational, diffusion, etc.
Definition at line 134 of file QMCFlags.h.
Referenced by QMCRun::branchWalkers(), QMCRun::calculatePopulationSizeBiasCorrectionFactor(), checkFlags(), QMCWalker::ID(), operator<<(), read_flags(), QMCRun::readXML(), QMCManager::readXML(), QMCWalker::reweight_walker(), QMCManager::run(), QMCRun::toXML(), QMCManager::writeEnergyResultsHeader(), and QMCManager::writeEnergyResultsSummary().
Whether to update one electron at a time (1) or all at once (0).
The discussions that I've found on the issue are pg 97, pgs 273-278 of the Hammond, Lester, Reynolds QMC book They say acceptance probability is "reduced considerably" if all electrons are moved at once.
pg 144 of the Umrigar "VMC basics and applications to atoms and molecules" from the NATO series, 1998... He also recommends moving one electron at a time.
On the other hand, if we move one electron per iteration, then the amount of system memory required scales with the number of walkers.
You will want to link LAPACK if you move all at once. If you move one at a time, we use a different method to update the inverse.
Definition at line 155 of file QMCFlags.h.
Referenced by checkFlags(), QMCManager::checkMaxStepsTerminationCriteria(), QMCWalker::createChildWalkers(), QMCWalkerData::initialize(), QMCWalker::initializePropagation(), operator<<(), QMCWalker::processPropagation(), read_flags(), QMCRun::step(), and QMCSlater::update_Ds().
Using QMCSurfer instead of QMCWalker.
This will let you "surf" the wavefunction.
Definition at line 162 of file QMCFlags.h.
Referenced by QMCWalker::processPropagation(), and read_flags().
Minimize the amount of output.
Definition at line 167 of file QMCFlags.h.
Referenced by QMCJastrowParameters::print(), qmcbeaver(), read_flags(), QMCManager::run(), and QMCManager::writeTimingData().
These are debugging flags.
They might be used in QMCSurfer.
Definition at line 172 of file QMCFlags.h.
Referenced by QMCSurfer::mainMenu(), and read_flags().
Definition at line 173 of file QMCFlags.h.
Referenced by QMCSurfer::mainMenu(), read_flags(), and QMCSurfer::scanEnergies().
string QMCFlags::temp_dir |
Scratch directory for the calculation.
Temporary files which will not be needed later will be written here.
Definition at line 179 of file QMCFlags.h.
Referenced by checkFlags(), QMCManager::initializeCalculationState(), operator<<(), read_flags(), set_filenames(), and QMCManager::writeCheckpoint().
Parallelization algorithm to use.
Definition at line 184 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
long QMCFlags::iseed |
Seed used in generating pseudo-random numbers for the calculation.
Definition at line 189 of file QMCFlags.h.
Referenced by QMCSCFJastrow::checkParameterDerivatives(), QMCManager::initialize(), operator<<(), and read_flags().
Implemented here is the method described in: S.
Chiesa, D.M. Ceperley, and S. Zhang, Phys. Rev. Lett. 94, 036404 (2005) This selects the method used in calculating the forces in a molecule. Choices are none, bare_hellmann_feynman, poly_hellmann_feynman, and slaterpoly_hellmann_feynman.
Definition at line 198 of file QMCFlags.h.
Referenced by QMCNuclearForces::calcCoefficients(), QMCWalker::calculateObservables(), checkFlags(), QMCSCFJastrow::evaluate(), QMCNuclearForces::evaluate(), QMCNuclearForces::getTemperTerm(), QMCWalkerData::initialize(), QMCWalker::initialize(), QMCSCFJastrow::initialize(), QMCRun::initialize(), QMCNuclearForces::initialize(), QMCManager::initialize(), operator<<(), QMCWalker::operator=(), qmcbeaver(), read_flags(), QMCWalker::resetFutureWalking(), and QMCManager::run().
string QMCFlags::sampling_method |
Diffusion Green's function to be used for the calculation.
Only importance sampled Green's functions should be used with DMC.
Definition at line 204 of file QMCFlags.h.
Referenced by QMCWalker::calculateReverseGreensFunction(), checkFlags(), QMCWalker::moveElectrons(), operator<<(), and read_flags().
Should the "quantum force" for the calculation be modified.
Modifications exist which, for example, can reduce the time-step bias in DMC.
Definition at line 211 of file QMCFlags.h.
Referenced by QMCSCFJastrow::calculate_Modified_Grad_PsiRatio(), operator<<(), and read_flags().
Should the "local energy" for the calculation be modified.
Modifications exist which, for example, can reduce the time-step bias in DMC.
Definition at line 218 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCWalker::reweight_walker().
If Umrigar's 1993 QMC algorithm (J.
Chem. Phys. 99, 2865(1993)) is used with the equal-electron modifications to the "local energy" and "quantum force," this is the parameter describing how the modifications will occur.
Definition at line 226 of file QMCFlags.h.
Referenced by QMCSCFJastrow::calculate_Modified_Grad_PsiRatio(), operator<<(), and read_flags().
Algorithm used to reweight the walkers after a time step.
Definition at line 231 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCWalker::reweight_walker().
double QMCFlags::trail_eps2 |
Trail's eps^2 value for his weighting method found in Alternative sampling for variational quantum Monte Carlo Phys Rev E 77, 016704 (2008).
He recommends that eps2 = variance
Definition at line 240 of file QMCFlags.h.
double QMCFlags::e_0 |
The e_0 value used in: Alternative sampling for variational quantum Monte Carlo Phys Rev E 77, 016704 (2008).
It is supposed to hold the best known VMC energy.
Definition at line 249 of file QMCFlags.h.
string QMCFlags::branching_method |
Method for branching walkers.
Possible methods include not branching, branching with all walkers getting unit weight, branching with walkers getting fractional weights, and similar methods.
Definition at line 257 of file QMCFlags.h.
Referenced by QMCRun::branchWalkers(), operator<<(), and read_flags().
If a walker's weight exceeds this threshold, it will branch when the branching method allows the walkers to have fractional weights.
Definition at line 264 of file QMCFlags.h.
Referenced by QMCWalker::branchRecommended(), QMCRun::nonunitWeightBranching(), operator<<(), read_flags(), and QMCWalker::reweight_walker().
double QMCFlags::fusion_threshold |
If a walker's weight is less than this threshold, it will fuse with another low weight walker when the branching method allows the walkers to have fractional weights.
Definition at line 271 of file QMCFlags.h.
Referenced by QMCRun::nonunitWeightBranching(), operator<<(), and read_flags().
1 if the global results are to be collected and broadcast to all nodes and used to calculate the estimated energy, trial energy, and effective time step.
Definition at line 278 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), QMCManager::readXML(), and QMCManager::run().
The interval at which the DMC ensemble is to be synchronized.
Definition at line 283 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
Parameter which can be used to prevent persistent configurations from interfering with a calculation.
A description of this is found in J. Chem. Phys. 99, 2865(1993).
Definition at line 290 of file QMCFlags.h.
Referenced by QMCWalker::acceptOrRejectMove(), QMCWalker::branchRecommended(), QMCWalker::calculateMoveAcceptanceProbability(), QMCWalker::calculateObservables(), QMCManager::equilibration_step(), operator<<(), and read_flags().
Warning verbosity level.
A higher integer should correspond to more warnings printed.
Definition at line 296 of file QMCFlags.h.
Referenced by QMCWalker::calculateObservables(), operator<<(), and read_flags().
double QMCFlags::rel_cutoff |
Controls.
Definition at line 301 of file QMCFlags.h.
Referenced by QMCWalker::calculateDerivatives(), QMCWalker::calculateObservables(), checkFlags(), operator<<(), QMCRun::randomlyInitializeWalkers(), read_flags(), and QMCWalker::reweight_walker().
Branching recommendations.
Definition at line 306 of file QMCFlags.h.
Referenced by QMCWalker::branchRecommended(), operator<<(), and read_flags().
double QMCFlags::dt |
Time step being used for the calculation.
Notes: There is sort of a balance here. If dt is large, then your error will be large. If dt is small, it will take longer to explore all the space.
Officially, you're supposed to estimate an extrapolation to dt = 0. That is, run it for (for example) 1e-3, then 1e-4, then 5e-5, etc until you can guess what the value at zero would be.
(AGA note, maybe allow this to be input as an array, and run all the calculations simultaneously or in succession...)
Definition at line 324 of file QMCFlags.h.
Referenced by QMCSCFJastrow::calculate_Modified_Grad_PsiRatio(), QMCRun::calculatePopulationSizeBiasCorrectionFactor(), QMCWalker::calculateReverseGreensFunction(), QMCWalker::calculateReverseGreensFunctionImportanceSampling(), QMCWalker::calculateReverseGreensFunctionUmrigar93ImportanceSampling(), checkFlags(), QMCManager::equilibration_step(), QMCWalkerData::getModifiedLocalEnergy(), QMCWalker::moveElectrons(), QMCWalker::moveElectronsImportanceSampling(), QMCWalker::moveElectronsNoImportanceSampling(), QMCWalker::moveElectronsUmrigar93ImportanceSampling(), operator<<(), read_flags(), QMCManager::run(), and QMCManager::updateEffectiveTimeStep().
double QMCFlags::accel_delta |
These two parameters are used by Umrigar's metropolis acceleration scheme, labeled here umrigar93_accelerated_sampling.
See QMCWalker for details.
Definition at line 331 of file QMCFlags.h.
Referenced by QMCWalker::calculateReverseGreensFunctionUmrigar93AcceleratedSampling(), checkFlags(), QMCWalker::moveElectronsUmrigar93AcceleratedSampling(), and read_flags().
double QMCFlags::accel_tm |
Definition at line 332 of file QMCFlags.h.
Referenced by QMCWalker::calculateReverseGreensFunctionUmrigar93AcceleratedSampling(), checkFlags(), QMCWalker::moveElectronsUmrigar93AcceleratedSampling(), and read_flags().
double QMCFlags::dt_run |
Time step to be used for the non-equilibration portion of the calculation.
Definition at line 338 of file QMCFlags.h.
Referenced by QMCManager::checkMaxTimeTerminationCriteria(), QMCManager::equilibration_step(), operator<<(), read_flags(), and QMCManager::run().
double QMCFlags::dt_effective |
Last updated value of the effective time step.
Definition at line 343 of file QMCFlags.h.
Referenced by QMCRun::calculatePopulationSizeBiasCorrectionFactor(), operator<<(), read_flags(), QMCWalker::reweight_walker(), QMCManager::updateEffectiveTimeStep(), QMCManager::updateTrialEnergy(), and QMCManager::writeEnergyResultsSummary().
Calculation terminates if the standard deviation of the energy falls below this threshold.
Notes: It's difficult to determine a priori how many time steps it will take to reach a given convergence, so it's often easier to control the length of the run by setting this to 0.0 and choosing an appropriate max_time_steps.
Definition at line 354 of file QMCFlags.h.
Referenced by QMCManager::checkConvergenceBasedTerminationCriteria(), operator<<(), and read_flags().
unsigned long QMCFlags::max_time_steps |
Maximum number of time steps the calculation will run for.
Set max_time_steps = 0 if you do not want it to be used as termination criteria.
Definition at line 361 of file QMCFlags.h.
Referenced by checkFlags(), QMCManager::checkMaxStepsTerminationCriteria(), operator<<(), QMCEigenSearch::optimize(), read_flags(), and QMCManager::run().
unsigned long QMCFlags::original_max_time_steps |
double QMCFlags::max_time |
Maximum amount of time the calculation will run for.
The number of steps that this will limit to will be max_time / dt. The main point of this parameter is to have one that will change automatically depending on dt. This makes it easier to use the same input files for different dt.
Set max_time = 0 if you do not want it to be used as termination criteria.
There is an ambiguitiy as far as what this parameter means when doing a parallel calculation. Should this limit the amount of time that each processor individually covers or the "collective" time that all processors cover? I've chosen the later since the former would involve a loss of parallelization efficiency.
Definition at line 380 of file QMCFlags.h.
Referenced by QMCManager::checkMaxTimeTerminationCriteria(), operator<<(), and read_flags().
This parameter is used in QMCRun to process several walkers simultaneously.
The value of this macro is how many walkers to treat at once. If this parameter is set to 1, then the code will perform the way it did previously. Also, the amount of memory the program requires is dependent on this because all the walkers have to be stored.
This used instead of the WALKERS_PER_PASS macro. DON"T CHANGE THIS IN THE MIDDLE OF A CALCULATION!
Definition at line 392 of file QMCFlags.h.
Referenced by QMCSlater::allocate(), QMCSlater::allocateIteration(), checkFlags(), getNumGPUWalkers(), QMCPotential_Energy::initialize(), QMCJastrow::initialize(), operator<<(), QMCRun::propagateWalkers(), read_flags(), and QMCSlater::~QMCSlater().
vector<int> QMCFlags::future_walking |
Implemented here is the future walking method described in: J.
Casulleras and J. Boronat, Phys. Rev. B 52, 3654 (1995)
This vector represents future walking projection time. Inspired by the decorrelation algorithm, it might be best to set this to something exponential like 1.0 2.0 4.0 8.0 16.0 32.0
A 0 is automatically included in the set since 0 represents no future walking. Preliminary results seem to indicate that future walking converges somewhere between 1 and 10 Hartrees^-1. For smaller dt, this requires more iterations.
In memory, we actually store it as block iteration length, since this is what QMCWalker uses.
Block Length = Time / dt
Definition at line 414 of file QMCFlags.h.
Referenced by QMCWalker::calculateObservables(), QMCWalker::initialize(), operator<<(), QMCPropertyArrays::QMCPropertyArrays(), read_flags(), QMCRun::readXML(), QMCRun::toXML(), and QMCPropertyArrays::toXML().
Definition at line 416 of file QMCFlags.h.
Referenced by checkFlags(), getNumGPUWalkers(), and read_flags().
Current number of walkers.
Definition at line 421 of file QMCFlags.h.
Referenced by checkFlags(), QMCRun::randomlyInitializeWalkers(), read_flags(), QMCRun::readXML(), and QMCManager::readXML().
Initial number of walkers.
Branching calculations will try to maintain this many walkers.
Definition at line 427 of file QMCFlags.h.
Referenced by QMCRun::ack_reconfiguration(), operator<<(), QMCRun::randomlyInitializeWalkers(), read_flags(), QMCManager::readXML(), QMCManager::run(), and QMCManager::writeTimingData().
Use interpolation of the radial basis functions in an attempt to speed up the calculation.
Definition at line 433 of file QMCFlags.h.
Referenced by QMCBasisFunction::initializeInterpolations(), operator<<(), and read_flags().
Number of grid points to use in interpolating the radial basis functions.
Definition at line 439 of file QMCFlags.h.
Referenced by QMCBasisFunction::initializeInterpolation(), operator<<(), and read_flags().
Smallest point to use in interpolating the radial basis functions.
Definition at line 444 of file QMCFlags.h.
Referenced by QMCBasisFunction::initializeInterpolation(), operator<<(), and read_flags().
This is the index of the Lebedev-Laikov grid to use for the run on each nucleus.
6 pts for psuedo_gridLevel = 0 14 pts for psuedo_gridLevel = 1 26 pts for psuedo_gridLevel = 2 38 pts for psuedo_gridLevel = 3 50 pts for psuedo_gridLevel = 4 etc see QMCMolecule::readPsuedoPotential for the other sizes
There are 32 levels to choose from but 14 or 26 points should be sufficient.
Definition at line 460 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), QMCInput::read(), and read_flags().
double QMCFlags::psuedo_cutoff |
If the local part of the psuedopotential is sufficiently small, then we can forgo the calculation of the nonlocal part.
Definition at line 466 of file QMCFlags.h.
Referenced by QMCPotential_Energy::evaluatePsuedoPotential(), operator<<(), and read_flags().
This is an internal flag, not an input flag.
It is turned on if we read any ECPs. It helps control memory allocation.
Definition at line 473 of file QMCFlags.h.
Referenced by QMCWalkerData::initialize(), QMCInput::read(), read_flags(), and QMCSlater::update_Ds().
Number of time steps taken between producing output.
Definition at line 478 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
1 if the accumulated statistics in the input checkpoint file are to be discarded, and 0 otherwise.
This is useful if you want to keep the (equilibrated) walker positions and weights, but don't want the accumulated statistics.
Definition at line 487 of file QMCFlags.h.
Referenced by QMCManager::initializeCalculationState(), operator<<(), qmcbeaver(), and read_flags().
unsigned long QMCFlags::equilibration_steps |
Number of time steps to take for the calculation to equilibrate.
Definition at line 492 of file QMCFlags.h.
Referenced by QMCWalker::calculateObservables(), checkFlags(), QMCManager::equilibration_step(), QMCWalker::initializePropagation(), QMCSurfer::mainMenu(), operator<<(), QMCManager::optimize(), QMCEigenSearch::optimize(), read_flags(), QMCWalker::reweight_walker(), QMCManager::run(), QMCRun::step(), and QMCManager::writeEnergyResultsSummary().
double QMCFlags::dt_equilibration |
Initial time step used during the equilibration of the calculation.
Notes: You might want to set this higher than you would set your regular dt since the point is to get the electrons in approximately the right place, not to make measurements.
Definition at line 501 of file QMCFlags.h.
Referenced by checkFlags(), QMCManager::equilibration_step(), operator<<(), read_flags(), and QMCManager::run().
Procedure used to modify the time step while the calculation is equilibrating.
Definition at line 507 of file QMCFlags.h.
Referenced by QMCManager::equilibration_step(), operator<<(), and read_flags().
unsigned int QMCFlags::CKAnnealingEquilibration1_parameter |
Number of time steps to use the CKAnnealingEquilibration1 before just equilibrating with the time step of the non-equilibration portion of the calculation.
The CKAnnealingEquilibration1 algorithm accepts essentially any proposed move. If importance sampling is used, this can rapidly reach a "good" section of configuration space.
Definition at line 517 of file QMCFlags.h.
Referenced by QMCManager::equilibration_step(), operator<<(), and read_flags().
1 if QMCEquilibrationArray is to be used.
This defines an array of QMCProperties objects, where the ith element starts collecting statistics on the (2^i)th iteration. The element with the lowest standard deviation for the total energy is chosen, automatically determining the correct number of equilibrtion steps. 0 if all the statistics are to be collected into one QMCProperties object.
Definition at line 528 of file QMCFlags.h.
Referenced by QMCRun::calculateObservables(), QMCManager::finalize(), QMCRun::getProperties(), QMCRun::initialize(), operator<<(), QMCManager::optimize(), read_flags(), QMCRun::readXML(), QMCManager::run(), QMCRun::toXML(), and QMCRun::zeroOut().
Number of time steps taken on the root processor before the results are collected from all other processors.
Notes: the longer this is, the longer it will take for the calculation to converge. The short this is, the more time the calculation will spend shuttling data back and forth.
This balance might be molecule specific since the communication time depends on how much data to send.
Definition at line 541 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
Number of time steps taken on the worker processors before checking for commands from the root processor.
Notes: if this is too large, then the root processor might end of spending a long time waiting. There isn't a significant penalty to making this small. Perhaps this is the sort of parameter we want to remain opaque to the user.
It should probably be pretty small relative to mpireduce_interval.
Definition at line 554 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
Number of time steps taken before this processor writes a checkpoint file.
Definition at line 560 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
1 if the calculation is to checkpoint, and 0 otherwise.
Definition at line 565 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), read_flags(), and QMCManager::run().
Use existing checkpoints when performing this calculation.
Definition at line 570 of file QMCFlags.h.
Referenced by QMCManager::initializeCalculationState(), operator<<(), and read_flags().
To make them smaller, write only the energy data to the checkpoint files.
Definition at line 575 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), QMCProperties::readXML(), and QMCProperties::toXML().
Number of time steps between printing the global properties as the calculation progresses.
Definition at line 581 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
1 if the global properties of the calculation are output as the calculation progresses, and 0 otherwise.
Definition at line 587 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
int QMCFlags::Natoms |
Number of atoms in the system.
Definition at line 592 of file QMCFlags.h.
Referenced by QMCBasisFunction::angularGrid(), QMCDansWalkerInitialization::assign_electrons_to_nuclei(), QMCBasisFunction::basisFunctionsOnGrid(), QMCWalker::calculateElectronDensities(), checkFlags(), QMCMikesJackedWalkerInitialization::electrons_and_radii(), QMCBasisFunction::evaluateBasisFunctions(), QMCMikesBetterWalkerInitialization::initializeBunchOfWalkersPosition(), QMCMikesJackedWalkerInitialization::initializeWalkerPosition(), QMCDansWalkerInitialization::initializeWalkerPosition(), QMCMikesBetterWalkerInitialization::ObjectiveFunctionForWalkers(), operator<<(), operator>>(), QMCInput::read(), and read_flags().
int QMCFlags::charge |
Charge of the system in atomic units.
(e.g. +1 if a neutral system looses an electron)
Definition at line 598 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), QMCInput::read(), and read_flags().
Number of orbitals in the trial wavefunction.
Definition at line 603 of file QMCFlags.h.
Referenced by operator<<(), QMCWavefunction::read(), QMCInput::read(), and read_flags().
Number of basis functions used to create the orbitals for the trial wavefunction.
Definition at line 609 of file QMCFlags.h.
Referenced by checkFlags(), QMCWalker::initialize(), operator<<(), QMCInput::read(), and read_flags().
Number of determinants in the SCF part of the wavefunction.
Definition at line 614 of file QMCFlags.h.
Referenced by operator<<(), QMCInput::read(), and read_flags().
Restricted or unrestricted trial functions can be used.
Definition at line 619 of file QMCFlags.h.
Referenced by checkFlags(), QMCWalker::initialize(), QMCRun::initializeFunction(), operator<<(), QMCInput::read(), and read_flags().
1 if the basis function density is to be calculated, 0 if not.
Definition at line 624 of file QMCFlags.h.
Referenced by QMCSlater::allocate(), QMCSCFJastrow::calculate_Psi_quantities(), QMCWalker::calculateObservables(), QMCSlater::evaluate(), QMCWalkerData::initialize(), QMCRun::initialize(), QMCManager::initialize(), operator<<(), QMCSlater::operator=(), qmcbeaver(), read_flags(), QMCManager::run(), and QMCSlater::~QMCSlater().
double QMCFlags::energy_trial |
Trial energy used for branching QMC calculations.
Definition at line 629 of file QMCFlags.h.
Referenced by QMCRun::calculatePopulationSizeBiasCorrectionFactor(), operator<<(), QMCEigenSearch::optimize(), read_flags(), QMCWalker::reweight_walker(), QMCRun::step(), QMCManager::updateTrialEnergy(), and QMCManager::writeEnergyResultsSummary().
double QMCFlags::energy_estimated |
Last estimate of the energy of the system.
Definition at line 634 of file QMCFlags.h.
Referenced by QMCObjectiveFunctionResult::calculate_umrigar88(), QMCObjectiveFunction::evaluate(), read_flags(), QMCWalker::reweight_walker(), QMCRun::step(), QMCManager::updateEstimatedEnergy(), and QMCManager::updateTrialEnergy().
Original estimate of the energy of the system read from the input file.
Definition at line 640 of file QMCFlags.h.
Referenced by QMCWalker::branchRecommended(), QMCWalker::calculateDerivatives(), QMCWalker::calculateObservables(), QMCRun::calculatePopulationSizeBiasCorrectionFactor(), QMCWalkerData::getModifiedLocalEnergy(), QMCRun::randomlyInitializeWalkers(), read_flags(), QMCWalker::reweight_walker(), and QMCRun::step().
1 if the correction for the population size bias is used, and 0 otherwise.
See J. Chem. Phys. 99, 2865(1993).
Definition at line 646 of file QMCFlags.h.
Referenced by QMCRun::calculatePopulationSizeBiasCorrectionFactor(), operator<<(), and read_flags().
1 if walkers from the calculation are written out to file, and 0 otherwise.
Notes: turning this on is required for wavefunction optimization. However, it can *easily* take up all of your compute time since it obviously requires harddrive access. If you select optimization, this will automatically be set to 1.
Also, depending on print_config_frequency and molecule size, these files can easily get up to be several gigabytes in size. The files can be output in ASCII or binary -- that option is set in QMCConfigIO.
Eventually, we'll probably add HDF5 functionality, but last I (AGA) checked they were still working on their C++ API.
Definition at line 664 of file QMCFlags.h.
Referenced by checkFlags(), QMCInput::openConfigFile(), operator<<(), read_flags(), and QMCManager::run().
Number of time steps between writing walkers out to file.
Definition at line 669 of file QMCFlags.h.
Referenced by operator<<(), QMCManager::optimize(), read_flags(), and QMCManager::run().
1 if the wavefunction is to be optimized during a VMC wavefunction, and 0 otherwise.
Currently, we can optimize: 1) EE Jastrow parameters 2) NE Jastrow parameters 3) CI coefficients
Notes: This will rerun the entire VMC calculation each time it optimizes. That is, your calculation time might go up by a factor of max_optimize_Psi_steps, depending on convergence.
Definition at line 684 of file QMCFlags.h.
Referenced by checkFlags(), QMCManager::initializeCalculationState(), operator<<(), QMCManager::optimize(), qmcbeaver(), QMCWavefunction::read(), read_flags(), QMCManager::run(), and QMCManager::writeRestart().
Number of vmc-optimize iterations used when optimizing the wavefunction.
Notes: Wavefunction convergence can take a long time, depending on your Jastrow parameters.
Definition at line 693 of file QMCFlags.h.
Referenced by operator<<(), qmcbeaver(), and read_flags().
Objective function parameter used in forming a barrier around a valid region of parameter space.
Not all objective functions use such a barrier. A valid region of parameter space is one where the wavefunction is not majorly different from the one used to generate the walkers when correlated sampling is used.
Definition at line 702 of file QMCFlags.h.
Referenced by QMCObjectiveFunctionResult::mikes_score_function(), operator<<(), and read_flags().
This is an internal flag.
Some flags will require us to calculate the analytic wavefunction parameter derivatives, others won't. QMCFlags will figure it out and set this parameter accordingly.
Definition at line 710 of file QMCFlags.h.
Referenced by QMCSCFJastrow::calculate_Psi_quantities(), QMCWalker::calculateDerivatives(), QMCWalker::calculateObservables(), checkFlags(), QMCSCFJastrow::checkParameterDerivatives(), QMCThreeBodyJastrow::collectForPair(), QMCJastrowElectronElectron::collectForPair(), QMCJastrowElectronNuclear::evaluate(), QMCJastrow::evaluate(), CambridgeThreeBodyCorrelationFunction::evaluate(), Cambridge2CorrelationFunction::evaluate(), CambridgeThreeBodyCorrelationFunction::initializeParameters(), QMCEigenSearch::optimize(), QMCCorrelatedSamplingVMCOptimization::optimize(), QMCThreeBodyJastrow::packageDerivatives(), read_flags(), QMCManager::run(), CambridgeThreeBodyCorrelationFunction::setElectron(), and QMCPropertyArrays::zeroOut().
These next flags control which part of the wavefunction we're optimizing.
If the optimization method is set to automatic, then these parameters will automatically be changed.
Definition at line 717 of file QMCFlags.h.
Referenced by QMCJastrowElectronElectron::collectForPair(), QMCJastrowParameters::getNumberEdnEdnParameters(), QMCJastrowParameters::getNumberEEParameters(), QMCJastrowParameters::getNumberEupEdnParameters(), QMCJastrowParameters::getNumberEupEupParameters(), operator<<(), and read_flags().
Definition at line 718 of file QMCFlags.h.
Referenced by QMCJastrowElectronNuclear::evaluate(), QMCJastrowParameters::getNumberNEParameters(), operator<<(), and read_flags().
Definition at line 719 of file QMCFlags.h.
Referenced by QMCThreeBodyJastrow::collectForPair(), QMCJastrowParameters::getJWParameters(), QMCJastrowParameters::getNumberNEdnEdnParameters(), QMCJastrowParameters::getNumberNEupEdnParameters(), QMCJastrowParameters::getNumberNEupEupParameters(), operator<<(), QMCThreeBodyJastrow::packageDerivatives(), read_flags(), and QMCJastrowParameters::setJWParameters().
Definition at line 720 of file QMCFlags.h.
Referenced by QMCSCFJastrow::calculate_CorrelatedSampling(), QMCWavefunction::getNumberCIParameters(), operator<<(), and read_flags().
Definition at line 721 of file QMCFlags.h.
Referenced by QMCSlater::allocate(), QMCSlater::allocateIteration(), QMCSCFJastrow::calculate_CorrelatedSampling(), QMCWavefunction::getNumberORParameters(), QMCWavefunction::getORParameters(), operator<<(), QMCWavefunction::read(), read_flags(), QMCWavefunction::setORParameters(), QMCSlater::update_Ds(), and QMCSlater::~QMCSlater().
Some Jastrows (e.g.
the cambridge ones) use cutoffs which can be difficult to optimize. this one parameter will let you decide whether to optimize them.
Definition at line 727 of file QMCFlags.h.
Referenced by CambridgeThreeBodyCorrelationFunction::evaluate(), Williamson2CorrelationFunction::get_p2_xa(), CambridgeThreeBodyCorrelationFunction::get_p2_xa(), Anderson2CorrelationFunction::get_p2_xa(), Williamson2CorrelationFunction::get_p3_xxa(), CambridgeThreeBodyCorrelationFunction::get_p3_xxa(), Anderson2CorrelationFunction::get_p3_xxa(), Williamson2CorrelationFunction::get_p_a(), CambridgeThreeBodyCorrelationFunction::get_p_a(), Anderson2CorrelationFunction::get_p_a(), Cambridge2CorrelationFunction::initializeParameters(), operator<<(), QMCEigenSearch::optimize(), CambridgeThreeBodyCorrelationFunction::print(), and read_flags().
Definition at line 729 of file QMCFlags.h.
Referenced by QMCJastrow::evaluate(), QMCPotential_Energy::evaluatePsuedoPotential(), QMCWalkerData::initialize(), QMCJastrow::initialize(), operator<<(), QMCJastrowParameters::operator=(), QMCJastrow::operator=(), QMCJastrowParameters::read(), read_flags(), and QMCJastrowParameters::setJWParameters().
If this flag is true, we allow the terms in the three body Jastrow that overlap the NE Jastrow to be nonzero.
Definition at line 735 of file QMCFlags.h.
Referenced by QMCThreeBodyCorrelationFunctionParameters::gaussianParamDepMatrix(), QMCThreeBodyCorrelationFunctionParameters::makeParamDepMatrix(), operator<<(), and read_flags().
If this flag is true, we allow the terms in the three body Jastrow that overlap the EE Jastrow to be nonzero.
Definition at line 741 of file QMCFlags.h.
Referenced by QMCThreeBodyCorrelationFunctionParameters::gaussianParamDepMatrix(), QMCThreeBodyCorrelationFunctionParameters::makeParamDepMatrix(), operator<<(), and read_flags().
double QMCFlags::a_diag |
This controls the initial value used by QMCLineSearch for modifying the hessian.
If this is large, then the optimization will be more like steepest descent, which is slower but more stable.
If it's smaller, then optimization should happen faster.
Definition at line 751 of file QMCFlags.h.
Referenced by QMCEigenSearch::get_a_diag(), CambridgeThreeBodyCorrelationFunction::initializeParameters(), Cambridge2CorrelationFunction::initializeParameters(), operator<<(), QMCLineSearch::optimize(), QMCEigenSearch::optimize(), QMCCorrelatedSamplingVMCOptimization::optimize(), qmcbeaver(), and read_flags().
double QMCFlags::ksi |
A parameter used by the Linearize step length method, must be between 0 and 1.
The method is most stable when ksi = 0.
Definition at line 759 of file QMCFlags.h.
Referenced by checkFlags(), QMCEigenSearch::getParameters(), operator<<(), and read_flags().
Objective function to optimize when optimizing a VMC wavefunction.
Definition at line 764 of file QMCFlags.h.
Referenced by QMCWalker::calculateDerivatives(), checkFlags(), QMCDerivativeProperties::getParameterGradient(), QMCDerivativeProperties::getParameterHamiltonian(), QMCDerivativeProperties::getParameterHessian(), QMCDerivativeProperties::getParameterOverlap(), QMCDerivativeProperties::getParameterValue(), operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), read_flags(), QMCObjectiveFunctionResult::set_score(), and QMCPropertyArrays::zeroOut().
Numerical optimization algorithm to use when optimizing a VMC wavefunction.
Notes: for QMC, if your wavefunction is an eigenfunction, your variance will be exactly zero. That's why you might want to minimize your variance.
If you set this to "automatic", then the code will attempt to select the parameter choices that probably work best. However, you might want to tweek some stuff anyway. For example, the a_diag parameter in QMCLineSearch.cpp the run lengths in QMcBeaver.cpp
Definition at line 783 of file QMCFlags.h.
Referenced by checkFlags(), QMCEigenSearch::get_a_diag(), operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), QMCLineSearch::optimize(), and read_flags().
Objective function to use when calculating derivatives for optimizing a VMC wavefunction.
Definition at line 789 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCObjectiveFunctionResult::set_score_for_derivative().
Maximum number of iterations to use when optimizing a VMC wavefunction with one set of correlated-sampling configurations.
Notes: This is how many iterations the optimization might take to converge a set of parameters over your wavefunction during a given optimization cycle.
Definition at line 799 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), QMCCorrelatedSamplingVMCOptimization::optimize(), and read_flags().
Tolerance used to determine when a numerical optimization with one set of correlated-sampling configurations has converged.
Definition at line 805 of file QMCFlags.h.
Referenced by operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), and read_flags().
Population size used when the CKGeneticAlgorithm1 is used to optimize a VMC wavefunction.
Definition at line 812 of file QMCFlags.h.
Referenced by operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), and read_flags().
Mutation rate used when the CKGeneticAlgorithm1 is used to optimize a VMC wavefunction.
Definition at line 819 of file QMCFlags.h.
Referenced by operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), and read_flags().
Standard deviation used when initializing the CKGeneticAlgorithm1.
Definition at line 825 of file QMCFlags.h.
Referenced by operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), and read_flags().
Algorithm used to determine the step length used when a line search algorithm is used to optimize a VMC wavefunction.
Definition at line 831 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), QMCOptimizationFactory::optimizationAlgorithmFactory(), and read_flags().
Parameter used for the logarithmic barrier when D.R.
Kent IV's algorithm for optimizing possibly singular Jastrow functions is used. See David Randall Kent IV's thesis, New Quantum Monte Carlo Algorithms to Efficiently Utilize Massively Parallel Computers, from the California Institute of Technology for more details. This parameter should be a small positive number.
Definition at line 841 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), QMCObjectiveFunctionResult::set_score(), and QMCObjectiveFunctionResult::set_score_for_derivative().
1 if up and down electrons use the same Jastrow parameters for their electron-nucleus terms, and 0 otherwise.
Definition at line 847 of file QMCFlags.h.
Referenced by QMCJastrowElectronNuclear::evaluate(), operator<<(), QMCJastrowParameters::read(), QMCInput::read(), read_flags(), QMCJastrowElectronElectron::updateAll(), and QMCJastrowElectronElectron::updateOne().
Definition at line 848 of file QMCFlags.h.
Referenced by QMCJastrowParameters::getJWParameters(), QMCJastrowParameters::getNumberNEdnEdnParameters(), QMCJastrowParameters::getNumberNEupEupParameters(), operator<<(), QMCThreeBodyJastrow::packageDerivatives(), QMCJastrowParameters::print(), QMCJastrowParameters::read(), read_flags(), and QMCJastrowParameters::setJWParameters().
1 if up and down electrons use the same orbital parameters.
Just because the input file used the same orbitals from the ckmf file doesn't mean they can't be optimized independently.
The restart file will have more orbitals than the input file.
Definition at line 858 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), QMCWavefunction::read(), and read_flags().
Just because the determinants shared orbitals when they were defined using the input file, that doesn't mean we can't optimize them independently.
The restart file will have more orbitals than the input file.
Definition at line 868 of file QMCFlags.h.
Referenced by operator<<(), QMCWavefunction::read(), and read_flags().
1 if Gaussian orbitals are to be replaced by exponentials that satisfy the electron nucleus cusp consitions near nuclei and 0 otherwise.
Definition at line 874 of file QMCFlags.h.
Referenced by checkFlags(), QMCSlater::evaluate(), QMCPotential_Energy::evaluatePsuedoPotential(), QMCSlater::initialize(), QMCPotential_Energy::initialize(), QMCSurfer::mainMenu(), operator<<(), QMCSlater::operator=(), QMCInput::read(), and read_flags().
string QMCFlags::energy_cutoff_type |
1 if we'll use the energy cutoff recommended by DePasquale, Rothstein, Vrbik JCP, 89, 3629, 1988 0 otherwise
Definition at line 881 of file QMCFlags.h.
Referenced by QMCWalkerData::getModifiedLocalEnergy(), operator<<(), and read_flags().
1 if every VMC calculation used to generated correlated-sampling configurations for an optimization equilibrates for the specified equilibration time, and 0 otherwise.
Definition at line 888 of file QMCFlags.h.
Referenced by operator<<(), qmcbeaver(), and read_flags().
1 if the first VMC calculation used to generated correlated-sampling configurations for an optimization equilibrates for the specified equilibration time, and 0 otherwise.
If a non-VMC calculation is performed, 1 if the calculation equilibrates, and 0 otherwise.
Definition at line 896 of file QMCFlags.h.
Referenced by QMCManager::initializeCalculationState(), operator<<(), read_flags(), and QMCManager::run().
Parameter used in controlling the number of walkers or total weight of all walkers during a branching QMC calculation.
Definition at line 902 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::updateTrialEnergy().
1 if the local energy for every walker for every time step is to be written out, and 0 otherwise.
Definition at line 908 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::run().
1 if the distance between each pair of electrons for every time step is to be collected and written out in a histogram, and 0 otherwise.
Definition at line 914 of file QMCFlags.h.
Referenced by QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::run().
double QMCFlags::max_pair_distance |
The maximum distance, in atomic units, of the histograms for the electron densities.
Definition at line 920 of file QMCFlags.h.
Referenced by QMCRun::initialize(), operator<<(), and read_flags().
1 if the x coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 927 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the x coordinate of the 2D parallel spin correlation diagram.
Definition at line 933 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the x coordinate of the 2D parallel spin correlation diagram.
Definition at line 939 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
1 if the y coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 946 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the y coordinate of the 2D parallel spin correlation diagram.
Definition at line 952 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the y coordinate of the 2D parallel spin correlation diagram.
Definition at line 958 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
1 if the z coordinates of each parallel spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 965 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the z coordinate of the 2D parallel spin correlation diagram.
Definition at line 971 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the z coordinate of the 2D parallel spin correlation diagram.
Definition at line 977 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
1 if the x coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 984 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the x coordinate of the 2D opposite spin correlation diagram.
Definition at line 990 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the x coordinate of the 2D opposite spin correlation diagram.
Definition at line 996 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
1 if the y coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 1003 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the y coordinate of the 2D opposite spin correlation diagram.
Definition at line 1009 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the y coordinate of the 2D opposite spin correlation diagram.
Definition at line 1015 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
1 if the z coordinates of each opposite spin pair of electrons are to be collected in a 2D histogram to make a correlation diagram, and 0 otherwise.
Definition at line 1022 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCManager::gatherHistograms(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The minimum value for the z coordinate of the 2D opposite spin correlation diagram.
Definition at line 1028 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
The maximum value for the z coordinate of the 2D opposite spin correlation diagram.
Definition at line 1034 of file QMCFlags.h.
Referenced by QMCRun::calculateElectronDensities(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::writeElectronDensityHistograms().
Are Chip and Mike cool? Answer: Yea Baby!
Definition at line 1039 of file QMCFlags.h.
Referenced by checkFlags(), operator<<(), and read_flags().
MPI rank of this processor.
Definition at line 1044 of file QMCFlags.h.
Referenced by QMCManager::finalizeOutputs(), QMCManager::gatherHistograms(), QMCWalker::ID(), QMCManager::initializeCalculationState(), QMCManager::initializeMPI(), QMCManager::initializeOutputs(), QMCManager::optimize(), QMCCorrelatedSamplingVMCOptimization::optimize(), qmcbeaver(), QMCFlags(), QMCManager::run(), QMCManager::sendAllProcessorsInputFileName(), set_filenames(), QMCInput::setMPIParameters(), QMCRun::step(), QMCManager::writeCheckpoint(), and QMCManager::writeRestart().
int QMCFlags::nprocs |
Number of processors used in this calculation.
Definition at line 1049 of file QMCFlags.h.
Referenced by QMCManager::initializeMPI(), QMCFlags(), QMCManager::sendAllProcessorsACommand(), and QMCInput::setMPIParameters().
Hartree-Fock calculations using QMC, aimed at producing basis independent reference wavefunctions.
Computes V_eff for each electron by averaging over multiple walkers and configurations.
Definition at line 1056 of file QMCFlags.h.
Referenced by QMCPotential_Energy::calc_P_ee(), QMCRun::initialize(), operator<<(), read_flags(), and QMCManager::run().
Number of electron samples to average v_eff over.
Definition at line 1061 of file QMCFlags.h.
Referenced by QMCHartreeFock::Initialize(), operator<<(), and read_flags().
Option to keep trial energy fixed to its original value.
Definition at line 1066 of file QMCFlags.h.
Referenced by operator<<(), read_flags(), and QMCManager::updateTrialEnergy().