pub struct Config {Show 63 fields
pub thresholds: Thresholds,
pub max_steps: usize,
pub max_step_size: f64,
pub reduced_factor: f64,
pub nprocs: usize,
pub mem: String,
pub charge: i32,
pub charge_b: Option<i32>,
pub mult_state_a: usize,
pub mult_state_b: usize,
pub method: String,
pub program: QMProgram,
pub run_mode: RunMode,
pub td_state_a: String,
pub td_state_b: String,
pub mp2: bool,
pub fix_de: f64,
pub scans: Vec<ScanSpec>,
pub bagel_model: String,
pub program_commands: HashMap<String, String>,
pub is_oniom: bool,
pub oniom_layer_info: Vec<String>,
pub charge_and_mult_oniom1: String,
pub charge_and_mult_oniom2: String,
pub basis_set: String,
pub solvent: String,
pub dispersion: String,
pub restart: bool,
pub custom_interface_file: String,
pub state_a: usize,
pub state_b: usize,
pub drive_coordinate: String,
pub drive_start: f64,
pub drive_end: f64,
pub drive_steps: usize,
pub drive_type: String,
pub drive_atoms: Vec<usize>,
pub use_gediis: bool,
pub use_hybrid_gediis: bool,
pub gediis_switch_rms: f64,
pub gediis_switch_step: f64,
pub switch_step: usize,
pub bfgs_rho: f64,
pub max_history: usize,
pub print_checkpoint: bool,
pub smart_history: bool,
pub print_level: usize,
pub use_robust_diis: bool,
pub gediis_variant: String,
pub gediis_blend_mode: String,
pub gdiis_cosine_check: String,
pub gdiis_coeff_check: String,
pub n_neg: usize,
pub gediis_sim_switch: f64,
pub hessian_method: HessianMethod,
pub step_reduction_multiplier: f64,
pub trust_reduction_factor: f64,
pub trust_increase_factor: f64,
pub trust_inc_threshold: f64,
pub trust_dec_threshold: f64,
pub trust_min_radius: f64,
pub trust_max_radius: f64,
pub steepest_descent_step: f64,
}Expand description
Complete configuration for an OpenMECP calculation.
The Config struct contains all parameters needed to run a MECP optimization,
including quantum chemistry settings, optimization parameters, constraints,
and advanced features. It can be parsed from an input file or created
programmatically.
§Required Fields
At minimum, you must specify:
method: Quantum chemistry method and basis set (e.g., “B3LYP/6-31G*”)program: QM program to use (Gaussian, ORCA, etc.)nprocs: Number of processorsmem: Memory allocation (e.g., “4GB”)charge: Molecular charge (shared for both states)mult_state_a,mult_state_b: Spin multiplicities for states A and B
§Optional Fields
Many optional parameters have sensible defaults:
max_steps: 100 optimization stepsmax_step_size: 0.1 Åthresholds: Standard convergence criteriause_gediis: false (uses GDIIS by default)
§Examples
use omecp::config::{Config, QMProgram, RunMode};
// Create a basic configuration
let mut config = Config::default();
config.method = "B3LYP/6-31G*".to_string();
config.program = QMProgram::Gaussian;
config.nprocs = 4;
config.mem = "4GB".to_string();
config.charge = 0;
config.mult_state_a = 1; // Singlet for state A
config.mult_state_b = 3; // Triplet for state BFields§
§thresholds: ThresholdsConvergence thresholds for optimization
max_steps: usizeMaximum number of optimization steps
max_step_size: f64Maximum step size in Angstrom (Å).
reduced_factor: f64Step size reduction factor for line search
nprocs: usizeNumber of processors for QM calculations
mem: StringMemory allocation (e.g., “4GB”, “8GB”)
charge: i32Charge for system (shared for both states)
charge_b: Option<i32>Separate charge for state B (None = same as charge)
mult_state_a: usizeSpin multiplicity for state A (2S+1, where S is total spin)
mult_state_b: usizeSpin multiplicity for state B
method: StringQuantum chemistry method and basis set (e.g., “B3LYP/6-31G*”)
program: QMProgramQuantum chemistry program to use
run_mode: RunModeExecution mode for the calculation
td_state_a: StringTD-DFT keywords for state A (Gaussian format)
td_state_b: StringTD-DFT keywords for state B (Gaussian format)
mp2: boolUse MP2 instead of DFT (if supported by QM program)
fix_de: f64Target energy difference in eV (for FixDE mode)
scans: Vec<ScanSpec>List of PES scans to perform
bagel_model: StringBAGEL model specification (for BAGEL program)
program_commands: HashMap<String, String>Custom command mappings for QM programs
is_oniom: boolEnable ONIOM (QM/MM) calculation
oniom_layer_info: Vec<String>ONIOM layer information (e.g., “H,L” for High, Low)
charge_and_mult_oniom1: StringONIOM charge/multiplicity for state 1
charge_and_mult_oniom2: StringONIOM charge/multiplicity for state 2
basis_set: StringBasis set specification (for programs that separate method/basis)
solvent: StringSolvent model specification
dispersion: StringDispersion correction
restart: boolEnable restart mode (read checkpoint file)
custom_interface_file: StringPath to custom QM interface JSON configuration
state_a: usizeState index for TD-DFT state A (0 = ground state, 1+ = excited states)
state_b: usizeState index for TD-DFT state B
drive_coordinate: StringReaction coordinate for path following
drive_start: f64Starting value for coordinate driving
drive_end: f64Ending value for coordinate driving
drive_steps: usizeNumber of steps for coordinate driving
drive_type: StringType of coordinate for driving (bond, angle, dihedral)
drive_atoms: Vec<usize>Atom indices for coordinate driving
use_gediis: boolUse GEDIIS optimizer instead of GDIIS (faster for difficult cases)
use_hybrid_gediis: boolEnable smart hybrid GEDIIS with production-grade adaptive weighting.
When true, uses smart_hybrid_gediis_step which automatically adjusts the blend between GDIIS and GEDIIS based on:
- Energy trend analysis (uphill detection)
- Linear regression deviation (oscillation detection)
- Scale-invariant relative metrics
The algorithm is calibrated on 1000+ real optimizations and provides robust convergence across diverse chemical systems.
Default: true (recommended - significantly more robust than fixed 50/50)
gediis_switch_rms: f64RMS gradient threshold for phase 1→2 switch in sequential hybrid (default: 0.005 Ha/Å)
Li & Frisch JCTC 2006: “switches to GEDIIS when the root-mean-square force of the latest point is smaller than 10⁻² au” 10⁻² Ha/Bohr ≈ 0.005 Ha/Å
gediis_switch_step: f64RMS displacement threshold for phase 2→3 switch in sequential hybrid (default: 0.001 Å)
Li & Frisch JCTC 2006: “when the root-mean-square RFO step of the latest point is less than 2.5×10⁻³ au, GDIIS is used until convergence” 2.5×10⁻³ Bohr ≈ 0.0013 Å, rounded to 0.001 Å
switch_step: usizeStep number at which to switch from BFGS to DIIS optimizers (default: 3)
- 0: Use DIIS from step 1 (no BFGS)
-
= max_steps: Use BFGS only (no DIIS)
- Other values: Switch from BFGS to DIIS at specified step
bfgs_rho: f64Scaling factor for BFGS step
max_history: usizeMaximum number of history entries for DIIS optimizers (GDIIS/GEDIIS)
Controls how many previous iterations are retained for interpolation in DIIS-based optimization methods. Larger values can improve convergence but use more memory. Default: 5.
print_checkpoint: boolEnable or disable checkpoint JSON file generation
When enabled (default: true), checkpoint files are saved during optimization to allow restarting calculations. When disabled, no checkpoint files are created. Supports: true/false, yes/no, 1/0
smart_history: boolEnable smart history management for DIIS optimizers
When enabled, removes the WORST point (rather than oldest) when history is full. Uses intelligent scoring based on:
- Energy difference from degeneracy
- Gradient norm
- Geometric redundancy
- MECP-specific gap penalties
Benefits: 20-30% faster convergence in some cases Default: false (uses traditional FIFO history)
Set to true if you want to experiment with smart history management.
print_level: usizeVerbosity level for output:
0: quiet (only convergence info)1: normal (step details)2: verbose (debug messages including DIIS internals)
use_robust_diis: boolUse the new Fortran-ported DIIS implementations (more robust).
When enabled, uses GdiisOptimizer and GediisOptimizer classes
which include:
- SR1 inverse matrix updates for GDIIS
- Multiple GEDIIS variants (RFO, Energy, Simultaneous)
- Cosine and coefficient validation
- Energy rise tracking and adaptive variant selection
Default: false (for backward compatibility)
gediis_variant: StringGEDIIS variant selection (only used when use_robust_diis = true).
Options:
- “auto”: Automatically select based on RMS error and energy trend
- “rfo”: RFO-DIIS using quadratic step overlaps
- “energy”: Energy-DIIS using gradient-coordinate products
- “simultaneous”: Combines Energy-DIIS with quadratic terms
When use_robust_diis = false (standard mode):
- “sequential” (default): Standard GDIIS/GEDIIS using mean inverse Hessian
- “blend”: New GDIIS_blend using mean true Hessian inversion
Also settable via use_gediis = blend / use_gediis = sequential in input.
Default: “auto”
gediis_blend_mode: StringBlend mode for hybrid GEDIIS+GDIIS blend methods (default: “fixed_sequential”).
"fixed": 50/50 fixed blend of GDIIS and EDIIS geometries."fixed_sequential": 50/50 blend far, pure GDIIS near convergence (default)."gradient": Gradient-weighted blend (EDIIS when far, GDIIS near minimum)."sequential": Phased switching (GDIIS → weighted → GDIIS).
Only used when use_hybrid_gediis = true and use_gediis = "blend".
gdiis_cosine_check: StringCosine check mode for GDIIS step validation.
Controls how the GDIIS step direction is validated against the last error vector. Options:
- “none”: No cosine check
- “zero”: CosLim = 0.0
- “standard”: CosLim = 0.71 (recommended)
- “variable”: Variable limit based on number of vectors
- “strict”: CosLim = 0.866
Default: “standard”
gdiis_coeff_check: StringCoefficient check mode for GDIIS validation.
Controls validation of DIIS coefficients to prevent excessive extrapolation. Options:
- “none”: No coefficient check
- “regular”: Standard coefficient check
- “force_recent”: Force recent vectors to have larger weight
- “combined”: Regular + ForceRecent
Default: “regular”
n_neg: usizeNumber of negative Hessian eigenvalues for saddle point search.
- 0: Minimum search (default)
- 1: Transition state search
-
1: Higher-order saddle points (rare)
Affects GEDIIS variant selection and TS scaling. Default: 0
gediis_sim_switch: f64RMS error threshold for GEDIIS variant switching (SimSw in Fortran).
When RMS error > this threshold, Energy-DIIS is preferred. When RMS error <= this threshold, RFO-DIIS is used.
Default: 0.0025 (matching Fortran)
hessian_method: HessianMethodHessian update method selection.
Options:
DirectPsb: Direct Hessian + PSB update (default, recommended). SolvesB × dk = -gvia LU decomposition. Better for saddle-point-like MECP.InverseBfgs: Inverse H⁻¹ + BFGS (Fortran-style).Bofill: Weighted Powell/MS for saddle points. Direct Hessian.Powell: Symmetric rank-one (SR1) update. Direct Hessian.BfgsPowellMix: Adaptive BFGS/Powell blend. Direct Hessian.
Default: DirectPsb
Compatibility:
- All direct Hessian methods (DirectPsb, Bofill, Powell, BfgsPowellMix)
are compatible with
use_gediis = blend. InverseBfgsis incompatible withuse_gediis = blend.
step_reduction_multiplier: f64Multiplier for RMS gradient threshold in step reduction activation. Step reduction triggers when history norm < rms_g × this value. Default: 10.0
trust_reduction_factor: f64Trust radius reduction factor when energy increases. Default: 0.5
trust_increase_factor: f64Trust radius increase factor when energy decreases. Default: 1.2
trust_inc_threshold: f64Energy increase threshold (Ha) for trust radius reduction. Default: 0.0001
trust_dec_threshold: f64Energy decrease threshold (Ha) for trust radius increase. Default: 0.0001
trust_min_radius: f64Minimum trust radius in Angstrom. Default: 0.01
trust_max_radius: f64Maximum trust radius in Angstrom. Default: 1.0
steepest_descent_step: f64Steepest descent step size in Angstrom for fallback when DIIS/BFGS fails. Default: 0.01
Implementations§
Source§impl Config
impl Config
Sourcepub fn gediis_variant_is_blend(&self) -> bool
pub fn gediis_variant_is_blend(&self) -> bool
Returns true when the GEDIIS variant is set to "blend".
"blend" activates the new GDIIS_blend implementation
with blended mean Hessian error vectors and trust-region step control.
Sourcepub fn validate_hessian_optimizer_compatibility(&self) -> Result<(), String>
pub fn validate_hessian_optimizer_compatibility(&self) -> Result<(), String>
Validates compatibility between the Hessian method and optimizer variant.
The blend optimizer (use_gediis = blend) requires a direct Hessian method
because it solves H·e = F via LU decomposition. Inverse BFGS is incompatible.
§Returns
Ok(())if the combination is validErr(String)with a descriptive error message if invalid
Trait Implementations§
Source§impl<'de> Deserialize<'de> for Config
impl<'de> Deserialize<'de> for Config
Source§fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
Auto Trait Implementations§
impl Freeze for Config
impl RefUnwindSafe for Config
impl Send for Config
impl Sync for Config
impl Unpin for Config
impl UnsafeUnpin for Config
impl UnwindSafe for Config
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self from the equivalent element of its
superset. Read more§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self is actually part of its subset T (and can be converted to it).§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset but without any property checks. Always succeeds.§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self to the equivalent element of its superset.