Skip to main content

sci_form/eht/
params.rs

1//! EHT parameter tables: VSIP values and Slater exponents per element/orbital.
2//!
3//! Hoffmann parameters for valence orbitals of common elements.
4
5use serde::Serialize;
6
7/// Definition of a single atomic orbital for EHT.
8#[derive(Debug, Clone, Serialize)]
9pub struct OrbitalDef {
10    /// Principal quantum number n.
11    pub n: u8,
12    /// Angular momentum quantum number l (0=s, 1=p, 2=d).
13    pub l: u8,
14    /// Label, e.g. "2s", "2p".
15    pub label: &'static str,
16    /// Valence State Ionization Potential (eV), used as H_ii.
17    pub vsip: f64,
18    /// Slater exponent ζ (bohr⁻¹).
19    pub zeta: f64,
20}
21
22/// Full EHT parameter set for one element.
23#[derive(Debug, Clone, Serialize)]
24pub struct EhtParams {
25    /// Atomic number.
26    pub z: u8,
27    /// Element symbol.
28    pub symbol: &'static str,
29    /// Valence orbital definitions.
30    pub orbitals: &'static [OrbitalDef],
31}
32
33/// Wolfsberg-Helmholtz constant (dimensionless).
34pub const DEFAULT_K: f64 = 1.75;
35
36// ─── Parameter tables ────────────────────────────────────────────────────────
37
38static H_ORBITALS: [OrbitalDef; 1] = [OrbitalDef {
39    n: 1,
40    l: 0,
41    label: "1s",
42    vsip: -13.6,
43    zeta: 1.3,
44}];
45
46static C_ORBITALS: [OrbitalDef; 2] = [
47    OrbitalDef {
48        n: 2,
49        l: 0,
50        label: "2s",
51        vsip: -21.4,
52        zeta: 1.625,
53    },
54    OrbitalDef {
55        n: 2,
56        l: 1,
57        label: "2p",
58        vsip: -11.4,
59        zeta: 1.625,
60    },
61];
62
63static N_ORBITALS: [OrbitalDef; 2] = [
64    OrbitalDef {
65        n: 2,
66        l: 0,
67        label: "2s",
68        vsip: -26.0,
69        zeta: 1.950,
70    },
71    OrbitalDef {
72        n: 2,
73        l: 1,
74        label: "2p",
75        vsip: -13.4,
76        zeta: 1.950,
77    },
78];
79
80static O_ORBITALS: [OrbitalDef; 2] = [
81    OrbitalDef {
82        n: 2,
83        l: 0,
84        label: "2s",
85        vsip: -32.3,
86        zeta: 2.275,
87    },
88    OrbitalDef {
89        n: 2,
90        l: 1,
91        label: "2p",
92        vsip: -14.8,
93        zeta: 2.275,
94    },
95];
96
97static F_ORBITALS: [OrbitalDef; 2] = [
98    OrbitalDef {
99        n: 2,
100        l: 0,
101        label: "2s",
102        vsip: -40.0,
103        zeta: 2.425,
104    },
105    OrbitalDef {
106        n: 2,
107        l: 1,
108        label: "2p",
109        vsip: -18.1,
110        zeta: 2.425,
111    },
112];
113
114static CL_ORBITALS: [OrbitalDef; 2] = [
115    OrbitalDef {
116        n: 3,
117        l: 0,
118        label: "3s",
119        vsip: -26.3,
120        zeta: 2.183,
121    },
122    OrbitalDef {
123        n: 3,
124        l: 1,
125        label: "3p",
126        vsip: -14.2,
127        zeta: 1.733,
128    },
129];
130
131static BR_ORBITALS: [OrbitalDef; 2] = [
132    OrbitalDef {
133        n: 4,
134        l: 0,
135        label: "4s",
136        vsip: -22.07,
137        zeta: 2.588,
138    },
139    OrbitalDef {
140        n: 4,
141        l: 1,
142        label: "4p",
143        vsip: -13.1,
144        zeta: 2.131,
145    },
146];
147
148static I_ORBITALS: [OrbitalDef; 2] = [
149    OrbitalDef {
150        n: 5,
151        l: 0,
152        label: "5s",
153        vsip: -18.0,
154        zeta: 2.679,
155    },
156    OrbitalDef {
157        n: 5,
158        l: 1,
159        label: "5p",
160        vsip: -12.7,
161        zeta: 2.322,
162    },
163];
164
165static S_ORBITALS: [OrbitalDef; 2] = [
166    OrbitalDef {
167        n: 3,
168        l: 0,
169        label: "3s",
170        vsip: -20.0,
171        zeta: 1.817,
172    },
173    OrbitalDef {
174        n: 3,
175        l: 1,
176        label: "3p",
177        vsip: -11.0,
178        zeta: 1.817,
179    },
180];
181
182static P_ORBITALS: [OrbitalDef; 2] = [
183    OrbitalDef {
184        n: 3,
185        l: 0,
186        label: "3s",
187        vsip: -18.6,
188        zeta: 1.600,
189    },
190    OrbitalDef {
191        n: 3,
192        l: 1,
193        label: "3p",
194        vsip: -14.0,
195        zeta: 1.600,
196    },
197];
198
199static SI_ORBITALS: [OrbitalDef; 2] = [
200    OrbitalDef {
201        n: 3,
202        l: 0,
203        label: "3s",
204        vsip: -17.3,
205        zeta: 1.383,
206    },
207    OrbitalDef {
208        n: 3,
209        l: 1,
210        label: "3p",
211        vsip: -9.2,
212        zeta: 1.383,
213    },
214];
215
216static B_ORBITALS: [OrbitalDef; 2] = [
217    OrbitalDef {
218        n: 2,
219        l: 0,
220        label: "2s",
221        vsip: -15.2,
222        zeta: 1.300,
223    },
224    OrbitalDef {
225        n: 2,
226        l: 1,
227        label: "2p",
228        vsip: -8.5,
229        zeta: 1.300,
230    },
231];
232
233/// All supported element parameter sets.
234static ALL_PARAMS: &[EhtParams] = &[
235    EhtParams {
236        z: 1,
237        symbol: "H",
238        orbitals: &H_ORBITALS,
239    },
240    EhtParams {
241        z: 5,
242        symbol: "B",
243        orbitals: &B_ORBITALS,
244    },
245    EhtParams {
246        z: 6,
247        symbol: "C",
248        orbitals: &C_ORBITALS,
249    },
250    EhtParams {
251        z: 7,
252        symbol: "N",
253        orbitals: &N_ORBITALS,
254    },
255    EhtParams {
256        z: 8,
257        symbol: "O",
258        orbitals: &O_ORBITALS,
259    },
260    EhtParams {
261        z: 9,
262        symbol: "F",
263        orbitals: &F_ORBITALS,
264    },
265    EhtParams {
266        z: 14,
267        symbol: "Si",
268        orbitals: &SI_ORBITALS,
269    },
270    EhtParams {
271        z: 15,
272        symbol: "P",
273        orbitals: &P_ORBITALS,
274    },
275    EhtParams {
276        z: 16,
277        symbol: "S",
278        orbitals: &S_ORBITALS,
279    },
280    EhtParams {
281        z: 17,
282        symbol: "Cl",
283        orbitals: &CL_ORBITALS,
284    },
285    EhtParams {
286        z: 35,
287        symbol: "Br",
288        orbitals: &BR_ORBITALS,
289    },
290    EhtParams {
291        z: 53,
292        symbol: "I",
293        orbitals: &I_ORBITALS,
294    },
295];
296
297/// Look up EHT parameters by atomic number.
298pub fn get_params(z: u8) -> Option<&'static EhtParams> {
299    ALL_PARAMS.iter().find(|p| p.z == z)
300}
301
302/// Count the total number of valence basis functions for an element.
303/// s → 1, p → 3, d → 5.
304pub fn num_basis_functions(z: u8) -> usize {
305    match get_params(z) {
306        Some(p) => p
307            .orbitals
308            .iter()
309            .map(|o| match o.l {
310                0 => 1,
311                1 => 3,
312                2 => 5,
313                _ => 0,
314            })
315            .sum(),
316        None => 0,
317    }
318}
319
320#[cfg(test)]
321mod tests {
322    use super::*;
323
324    #[test]
325    fn test_hydrogen_params() {
326        let p = get_params(1).expect("H params");
327        assert_eq!(p.symbol, "H");
328        assert_eq!(p.orbitals.len(), 1);
329        assert!((p.orbitals[0].vsip - (-13.6)).abs() < 1e-10);
330        assert!((p.orbitals[0].zeta - 1.3).abs() < 1e-10);
331    }
332
333    #[test]
334    fn test_carbon_params() {
335        let p = get_params(6).expect("C params");
336        assert_eq!(p.symbol, "C");
337        assert_eq!(p.orbitals.len(), 2);
338        // 2s + 2px,2py,2pz = 4 basis functions
339        assert_eq!(num_basis_functions(6), 4);
340    }
341
342    #[test]
343    fn test_oxygen_params() {
344        let p = get_params(8).expect("O params");
345        assert_eq!(p.symbol, "O");
346        assert!((p.orbitals[0].vsip - (-32.3)).abs() < 1e-10);
347        assert!((p.orbitals[1].vsip - (-14.8)).abs() < 1e-10);
348    }
349
350    #[test]
351    fn test_all_elements_have_params() {
352        for z in [1, 5, 6, 7, 8, 9, 14, 15, 16, 17, 35, 53] {
353            assert!(get_params(z).is_some(), "Missing params for Z={}", z);
354        }
355    }
356
357    #[test]
358    fn test_unknown_element_returns_none() {
359        assert!(get_params(118).is_none());
360    }
361
362    #[test]
363    fn test_basis_function_count() {
364        assert_eq!(num_basis_functions(1), 1); // H: 1s
365        assert_eq!(num_basis_functions(6), 4); // C: 2s + 2p(3)
366        assert_eq!(num_basis_functions(8), 4); // O: 2s + 2p(3)
367    }
368}