| |
| """ |
| Compute AO brackets (U/V matrices + eigenvalues) for EPC from HPRO and DeepH-E3. |
| |
| Pipeline: |
| 1. HPRO reconstruction for each EPC structure (scf_0 + group_1..12) |
| 2. DeepH-E3 inference on each structure |
| 3. AO bracket computation: U[n,m](k) = c_n^{d,dag} · S^{d,0}(k) · c_m^0(k) |
| 4. Save ao_brackets_hpro.npz and ao_brackets_e3.npz in displacements/scf_0/ |
| |
| Run after diamond.jl (create + run steps), before diamond_ml.jl. |
| |
| Usage: |
| python ml_epc.py [--hpro] [--e3] [--brackets] (default: all steps) |
| python ml_epc.py --hpro # only HPRO reconstruction |
| python ml_epc.py --e3 # only DeepH-E3 inference |
| python ml_epc.py --brackets # only AO bracket computation |
| """ |
| import argparse |
| import glob |
| import json |
| import os |
| import subprocess |
| import sys |
|
|
| import h5py |
| import numpy as np |
| from ase.io import read as ase_read |
| from scipy import linalg as la |
|
|
| SCRIPT_DIR = os.path.dirname(os.path.abspath(__file__)) |
| DISP_DIR = os.path.join(SCRIPT_DIR, 'displacements') |
| DIAMOND_DIR = os.path.dirname(SCRIPT_DIR) |
| AOBASIS_DIR = os.path.join(DIAMOND_DIR, 'aobasis') |
| PSEUDOS_DIR = os.path.join(DIAMOND_DIR, 'pseudos') |
| PARAMS_PATH = os.path.join(DIAMOND_DIR, '1_data_prepare', 'params.json') |
| RESULTS_DIR = os.path.join(DIAMOND_DIR, '2_training', 'hamiltonian', 'results') |
|
|
| NATOMS = 2 |
| NBANDS = 8 |
| K_MESH = (6, 6, 6) |
|
|
| |
| _b = 3.567 / 2.0 |
| LAT_MAT = np.array([[0, _b, _b], [_b, 0, _b], [_b, _b, 0]]) |
| LAT_INV = np.linalg.inv(LAT_MAT) |
|
|
|
|
| |
| |
| |
|
|
| def generate_kgrid(n1, n2, n3): |
| kpts = [] |
| for i in range(n1): |
| for j in range(n2): |
| for k in range(n3): |
| kpts.append([i / n1, j / n2, k / n3]) |
| return np.array(kpts) |
|
|
|
|
| KGRID = generate_kgrid(*K_MESH) |
| NK = len(KGRID) |
|
|
|
|
| def read_h5_dict(fn): |
| """Read DeePH/HPRO HDF5 file → dict {(R0,R1,R2,ia,ja): ndarray}.""" |
| d = {} |
| with h5py.File(fn, 'r') as f: |
| for key in f.keys(): |
| nk = tuple(map(int, key[1:-1].split(','))) |
| d[nk] = np.array(f[key]) |
| return d |
|
|
|
|
| def get_orbital_info(aodir): |
| """Return (site_norbits, norbits, csum) from orbital_types.dat.""" |
| with open(os.path.join(aodir, 'orbital_types.dat')) as f: |
| site_norbits = np.array( |
| [sum(2 * l + 1 for l in map(int, line.split())) for line in f] |
| ) |
| norbits = int(sum(site_norbits)) |
| csum = np.cumsum(site_norbits) |
| return site_norbits, norbits, csum |
|
|
|
|
| def build_HR(h5d, s5d, site_norbits, norbits, csum): |
| """Assemble full H(R) and S(R) matrices from HDF5 dicts.""" |
| H_R, S_R = {}, {} |
| for key in h5d: |
| R, ai, aj = key[:3], key[3] - 1, key[4] - 1 |
| if R not in H_R: |
| H_R[R] = np.zeros((norbits, norbits), dtype=np.complex128) |
| S_R[R] = np.zeros((norbits, norbits), dtype=np.complex128) |
| i0 = csum[ai] - site_norbits[ai] |
| j0 = csum[aj] - site_norbits[aj] |
| ni, nj = site_norbits[ai], site_norbits[aj] |
| H_R[R][i0:i0 + ni, j0:j0 + nj] = h5d[key] |
| if key in s5d: |
| S_R[R][i0:i0 + ni, j0:j0 + nj] = s5d[key] |
| return H_R, S_R |
|
|
|
|
| def build_Hk_Sk(H_R, S_R, kf, norbits): |
| """Fourier-transform H(R), S(R) to k-space with Hermitianization.""" |
| Hk = np.zeros((norbits, norbits), dtype=np.complex128) |
| Sk = np.zeros((norbits, norbits), dtype=np.complex128) |
| for R in H_R: |
| ph = np.exp(2j * np.pi * np.dot(kf, R)) |
| Hk += H_R[R] * ph |
| if R in S_R: |
| Sk += S_R[R] * ph |
| return 0.5 * (Hk + Hk.conj().T), 0.5 * (Sk + Sk.conj().T) |
|
|
|
|
| def compute_eigenvectors(H_R, S_R, kpts, norbits, nbands=NBANDS): |
| """Solve generalized eigenproblem at each k-point, return (eigs, evecs) lists.""" |
| all_eigs, all_evecs = [], [] |
| for kf in kpts: |
| Hk, Sk = build_Hk_Sk(H_R, S_R, kf, norbits) |
| eigs, evecs = la.eigh(Hk, Sk) |
| idx = np.argsort(np.real(eigs)) |
| all_eigs.append(np.real(eigs[idx[:nbands]])) |
| all_evecs.append(evecs[:, idx[:nbands]]) |
| return all_eigs, all_evecs |
|
|
|
|
| def build_cross_overlap_R(S0_raw, Sd_raw, displaced_atom_1based, site_norbits, norbits, csum): |
| """Build S^{d,0}(R): use displaced-atom rows from Sd, others from S0.""" |
| S_cross = {} |
| all_keys = set(S0_raw.keys()) | set(Sd_raw.keys()) |
| for key in all_keys: |
| R = key[:3] |
| ai, aj = key[3] - 1, key[4] - 1 |
| if R not in S_cross: |
| S_cross[R] = np.zeros((norbits, norbits), dtype=np.complex128) |
| i0 = csum[ai] - site_norbits[ai] |
| j0 = csum[aj] - site_norbits[aj] |
| ni, nj = site_norbits[ai], site_norbits[aj] |
| if key[3] == displaced_atom_1based: |
| if key in Sd_raw: |
| S_cross[R][i0:i0 + ni, j0:j0 + nj] = Sd_raw[key] |
| else: |
| if key in S0_raw: |
| S_cross[R][i0:i0 + ni, j0:j0 + nj] = S0_raw[key] |
| return S_cross |
|
|
|
|
| _ORBPAIR_CACHE = None |
|
|
|
|
| def _get_orbpairs(): |
| """Initialize (and cache) HPRO OrbPair objects for Carbon DZP basis.""" |
| global _ORBPAIR_CACHE |
| if _ORBPAIR_CACHE is None: |
| from HPRO.orbutils import parse_siesta_ion, OrbPair, LinearRGD, grid_R2G |
| from HPRO.constants import AOFT_QGRID_DEN |
| ion_file = os.path.join(AOBASIS_DIR, 'C.ion') |
| norb_rad, phirgrids = parse_siesta_ion(ion_file) |
| Ecut = 50 |
| grid_nq = int(np.sqrt(Ecut) * AOFT_QGRID_DEN) |
| gridQ = LinearRGD(0, np.sqrt(2 * Ecut), grid_nq) |
| phiQlist = [grid_R2G(gridQ, g) for g in phirgrids] |
| orbpairs = [] |
| for jorb in range(norb_rad): |
| for iorb in range(norb_rad): |
| r1, r2 = phirgrids[iorb].rcut, phirgrids[jorb].rcut |
| orbpairs.append(OrbPair(phiQlist[iorb], phiQlist[jorb], r1 + r2, index=1)) |
| orbslices = [0] |
| for g in phirgrids: |
| orbslices.append(orbslices[-1] + 2 * g.l + 1) |
| |
| perm_full = np.array( |
| [0, 1] |
| + [2 + i for i in [2, 0, 1]] |
| + [5 + i for i in [2, 0, 1]] |
| + [8 + i for i in [2, 4, 0, 3, 1]] |
| ) |
| _ORBPAIR_CACHE = (orbpairs, orbslices, perm_full, norb_rad) |
| return _ORBPAIR_CACHE |
|
|
|
|
| def _compute_overlap_block(Rvec_cart_bohr, orbpairs, orbslices, norb_rad): |
| """Compute nao×nao overlap at given inter-atomic vector (Bohr) via OrbPair.""" |
| nao = orbslices[-1] |
| S = np.zeros((nao, nao)) |
| Rv = np.array(Rvec_cart_bohr).reshape(1, 3) |
| idx = 0 |
| for jorb in range(norb_rad): |
| for iorb in range(norb_rad): |
| s1 = slice(orbslices[iorb], orbslices[iorb + 1]) |
| s2 = slice(orbslices[jorb], orbslices[jorb + 1]) |
| S[s1, s2] = orbpairs[idx].calc(Rv)[0] |
| idx += 1 |
| return S |
|
|
|
|
| def _get_delta_cart_ang(scf0_dir, struct_dir_d, atom_1based): |
| """Cartesian displacement (Ang) of atom_1based between two QE scf.in files.""" |
| s0 = ase_read(os.path.join(scf0_dir, 'scf.in'), format='espresso-in') |
| sd = ase_read(os.path.join(struct_dir_d, 'scf.in'), format='espresso-in') |
| return sd.get_positions()[atom_1based - 1] - s0.get_positions()[atom_1based - 1] |
|
|
|
|
| def build_cross_overlap_R_exact(S0_raw, Sd_raw, displaced_atom_1based, |
| site_norbits, norbits, csum, delta_cart_ang): |
| """ |
| Exact cross-overlap S^{d,0}(R): same as build_cross_overlap_R but replaces |
| the (kappa,kappa) blocks with OrbPair evaluated at R_cart - delta_cart. |
| |
| For the self-overlap block of the displaced atom: |
| S^{d,0}_{kk}(R) = <phi(r - tau_k - delta) | phi(r - tau_k - R*lat)> |
| = OrbPair(R_cart - delta) |
| instead of OrbPair(R_cart) used in the approximate version. |
| """ |
| from HPRO.constants import bohr2ang |
| orbpairs, orbslices, perm_full, norb_rad = _get_orbpairs() |
| nao_per_atom = orbslices[-1] |
| P = np.eye(nao_per_atom)[perm_full] |
|
|
| kappa = displaced_atom_1based |
| ai_kappa = kappa - 1 |
| i0 = csum[ai_kappa] - site_norbits[ai_kappa] |
| ni = site_norbits[ai_kappa] |
| delta_bohr = np.array(delta_cart_ang) / bohr2ang |
|
|
| S_cross = build_cross_overlap_R(S0_raw, Sd_raw, kappa, site_norbits, norbits, csum) |
|
|
| for R in list(S_cross.keys()): |
| key_kk = R + (kappa, kappa) |
| if key_kk not in S0_raw and key_kk not in Sd_raw: |
| continue |
| R_cart_bohr = np.array(R) @ LAT_MAT / bohr2ang |
| R_cross_bohr = R_cart_bohr - delta_bohr |
| S_wiki = _compute_overlap_block(R_cross_bohr, orbpairs, orbslices, norb_rad) |
| S_cross[R][i0:i0 + ni, i0:i0 + ni] = P @ S_wiki @ P.T |
|
|
| return S_cross |
|
|
|
|
| def build_Sk_cross(S_cross_R, kf, norbits): |
| """Fourier-transform cross-overlap (no Hermitianization).""" |
| Sk = np.zeros((norbits, norbits), dtype=np.complex128) |
| for R in S_cross_R: |
| ph = np.exp(2j * np.pi * np.dot(kf, R)) |
| Sk += S_cross_R[R] * ph |
| return Sk |
|
|
|
|
| def compute_brackets(c_disp, c_undisp, S_cross_k, nbands=NBANDS): |
| """U[n,m] = c_n^{d,dag} · S^{d,0}(k) · c_m^0.""" |
| return c_disp[:, :nbands].conj().T @ S_cross_k @ c_undisp[:, :nbands] |
|
|
|
|
| |
| |
| |
|
|
| def run_hpro_for_structure(struct_dir, label): |
| """Run HPRO PW2AOkernel reconstruction for one EPC structure.""" |
| recon_dir = os.path.join(struct_dir, 'reconstruction') |
| aodir = os.path.join(recon_dir, 'aohamiltonian') |
| out_h5 = os.path.join(aodir, 'hamiltonians.h5') |
|
|
| if os.path.exists(out_h5): |
| print(f' [{label}] HPRO already done') |
| return True |
|
|
| |
| vscdir = os.path.join(struct_dir, 'tmp') |
| if not os.path.exists(os.path.join(vscdir, 'VSC')): |
| print(f' [{label}] WARNING: VSC not found at {vscdir}/VSC — run diamond.jl first') |
| return False |
|
|
| params = json.load(open(PARAMS_PATH)) |
| ecutwfn = params['hpro']['ecutwfn'] |
|
|
| os.makedirs(recon_dir, exist_ok=True) |
| orig_dir = os.getcwd() |
| os.chdir(recon_dir) |
| try: |
| from HPRO import PW2AOkernel |
| print(f' [{label}] Running HPRO reconstruction...') |
| kernel = PW2AOkernel( |
| lcao_interface='siesta', |
| lcaodata_root=AOBASIS_DIR, |
| hrdata_interface='qe-bgw', |
| vscdir=vscdir, |
| upfdir=PSEUDOS_DIR, |
| ecutwfn=ecutwfn, |
| ) |
| kernel.run_pw2ao_rs('./aohamiltonian') |
| print(f' [{label}] HPRO done → {aodir}') |
| finally: |
| os.chdir(orig_dir) |
| return True |
|
|
|
|
| def run_all_hpro(): |
| struct_labels = [('scf_0', os.path.join(DISP_DIR, 'scf_0'))] |
| for i in range(1, NATOMS * 6 + 1): |
| struct_labels.append((f'group_{i}', os.path.join(DISP_DIR, f'group_{i}'))) |
|
|
| print('=' * 60) |
| print('Step 1: HPRO reconstruction') |
| print('=' * 60) |
| for label, struct_dir in struct_labels: |
| if not os.path.isdir(struct_dir): |
| print(f' [{label}] not found — run diamond.jl (create + run) first') |
| continue |
| run_hpro_for_structure(struct_dir, label) |
|
|
|
|
| |
| |
| |
|
|
| def find_latest_model(): |
| dirs = sorted(glob.glob(os.path.join(RESULTS_DIR, '*'))) |
| dirs = [d for d in dirs if os.path.isdir(d) |
| and os.path.exists(os.path.join(d, 'best_model.pkl'))] |
| return dirs[-1] if dirs else None |
|
|
|
|
| def run_e3_inference(struct_dir, label, model_dir, params): |
| """Run DeepH-E3 inference for one structure, output hamiltonians_pred_e3.h5.""" |
| aodir = os.path.join(struct_dir, 'reconstruction', 'aohamiltonian') |
| pred_h5 = os.path.join(aodir, 'hamiltonians_pred_e3.h5') |
|
|
| if os.path.exists(pred_h5): |
| print(f' [{label}] DeepH-E3 prediction already exists') |
| return True |
| if not os.path.exists(os.path.join(aodir, 'hamiltonians.h5')): |
| print(f' [{label}] SKIP — no hamiltonians.h5 (run HPRO first)') |
| return False |
|
|
| infer_dir = os.path.join(struct_dir, 'infer_e3') |
| dataset_dir = os.path.join(infer_dir, 'dataset') |
| graph_dir = os.path.join(infer_dir, 'graph') |
| output_dir = os.path.join(infer_dir, 'output') |
| ini_path = os.path.join(infer_dir, 'eval.ini') |
|
|
| |
| dest = os.path.join(dataset_dir, '00') |
| os.makedirs(dest, exist_ok=True) |
| for fname in os.listdir(aodir): |
| link = os.path.join(dest, fname) |
| if not os.path.exists(link): |
| os.symlink(os.path.join(aodir, fname), link) |
| os.makedirs(graph_dir, exist_ok=True) |
| os.makedirs(output_dir, exist_ok=True) |
|
|
| |
| t = params.get('hamiltonian', {}) |
| conda_env = t.get('conda_env', 'deeph') |
| deeph_e3_dir = t.get('deeph_e3_dir', '/home/apolyukhin/Development/DeepH-E3') |
| device = t.get('device', 'cuda') |
|
|
| ini = f"""; DeepH-E3 eval config — generated by ml_epc.py |
| |
| [basic] |
| device = {device} |
| dtype = float |
| trained_model_dir = {model_dir} |
| output_dir = {output_dir} |
| target = hamiltonian |
| inference = False |
| test_only = False |
| |
| [data] |
| graph_dir = |
| DFT_data_dir = |
| processed_data_dir = {dataset_dir} |
| save_graph_dir = {graph_dir} |
| target_data = hamiltonian |
| dataset_name = diamond_qe_e3_epc_{label} |
| get_overlap = False |
| """ |
| with open(ini_path, 'w') as f: |
| f.write(ini) |
|
|
| |
| launcher = os.path.join(infer_dir, '_eval_launcher.py') |
| with open(launcher, 'w') as f: |
| f.write(f"""import sys, torch |
| torch.serialization.add_safe_globals([slice]) |
| try: |
| from torch_geometric.data.data import DataEdgeAttr, DataTensorAttr |
| from torch_geometric.data.storage import GlobalStorage |
| torch.serialization.add_safe_globals([DataEdgeAttr, DataTensorAttr, GlobalStorage]) |
| except ImportError: |
| pass |
| sys.path.insert(0, '{deeph_e3_dir}') |
| from deephe3 import DeepHE3Kernel |
| kernel = DeepHE3Kernel() |
| kernel.eval('{ini_path}') |
| """) |
|
|
| conda_base = params['paths']['conda_base'] |
| activate = (f'source {conda_base}/etc/profile.d/conda.sh' |
| f' && conda activate {conda_env}') |
| print(f' [{label}] Running DeepH-E3 inference...') |
| subprocess.run(['bash', '-c', f'{activate} && python {launcher}'], check=True) |
|
|
| |
| pred_out = os.path.join(output_dir, '00', 'hamiltonians_pred.h5') |
| if os.path.exists(pred_out): |
| import shutil |
| shutil.copy(pred_out, pred_h5) |
| print(f' [{label}] Inference done → {pred_h5}') |
| return True |
| else: |
| print(f' [{label}] WARNING: {pred_out} not found after inference') |
| return False |
|
|
|
|
| def run_all_e3(): |
| params = json.load(open(PARAMS_PATH)) |
| model_dir = find_latest_model() |
| if model_dir is None: |
| print('ERROR: No trained model found in 2_training/hamiltonian/results/') |
| sys.exit(1) |
| print(f'Using model: {os.path.basename(model_dir)}') |
|
|
| struct_labels = [('scf_0', os.path.join(DISP_DIR, 'scf_0'))] |
| for i in range(1, NATOMS * 6 + 1): |
| struct_labels.append((f'group_{i}', os.path.join(DISP_DIR, f'group_{i}'))) |
|
|
| print('=' * 60) |
| print('Step 2: DeepH-E3 inference') |
| print('=' * 60) |
| for label, struct_dir in struct_labels: |
| if not os.path.isdir(struct_dir): |
| continue |
| run_e3_inference(struct_dir, label, model_dir, params) |
|
|
|
|
| |
| |
| |
|
|
| def load_structure(struct_dir, ham_file): |
| """Load H_raw, S_raw, orbital info for one structure.""" |
| aodir = os.path.join(struct_dir, 'reconstruction', 'aohamiltonian') |
| H_raw = read_h5_dict(os.path.join(aodir, ham_file)) |
| S_raw = read_h5_dict(os.path.join(aodir, 'overlaps.h5')) |
| site_norbits, norbits, csum = get_orbital_info(aodir) |
| return H_raw, S_raw, site_norbits, norbits, csum |
|
|
|
|
| def compute_ao_brackets(source): |
| """ |
| Compute U/V brackets and eigenvalues for all 6 displacement directions. |
| |
| source: 'hpro' or 'e3' |
| Saves to displacements/scf_0/ao_brackets_{source}.npz |
| """ |
| ham_file = 'hamiltonians.h5' if source == 'hpro' else 'hamiltonians_pred_e3.h5' |
| label = source.upper() |
| out_file = os.path.join(DISP_DIR, 'scf_0', f'ao_brackets_{source}.npz') |
|
|
| if os.path.exists(out_file): |
| print(f' [{label}] ao_brackets_{source}.npz already exists') |
| return |
|
|
| print(f'\n{"=" * 60}') |
| print(f'Computing AO brackets ({label})') |
| print(f'{"=" * 60}') |
|
|
| |
| scf0_dir = os.path.join(DISP_DIR, 'scf_0') |
| H0_raw, S0_raw, site_norbits, norbits, csum = load_structure(scf0_dir, ham_file) |
| H0_R, S0_R = build_HR(H0_raw, S0_raw, site_norbits, norbits, csum) |
|
|
| print(f' Undisplaced: norbits={norbits}') |
| print(f' Computing eigenvectors at {NK} k-points...') |
| eigs_0, evecs_0 = compute_eigenvectors(H0_R, S0_R, KGRID, norbits) |
|
|
| n_disp_dirs = NATOMS * 3 |
| U_arr = np.zeros((n_disp_dirs, NK, NK, NBANDS, NBANDS), dtype=np.complex128) |
| V_arr = np.zeros((n_disp_dirs, NK, NK, NBANDS, NBANDS), dtype=np.complex128) |
| ek_arr = np.array(eigs_0) |
| ep_arr = np.zeros((n_disp_dirs, NK, NBANDS)) |
| epm_arr = np.zeros((n_disp_dirs, NK, NBANDS)) |
|
|
| for d in range(n_disp_dirs): |
| gp = 2 * d + 1 |
| gm = 2 * d + 2 |
| atom = d // 3 + 1 |
| print(f'\n d={d} (atom {atom}, groups {gp}/{gm})...') |
|
|
| |
| gp_dir = os.path.join(DISP_DIR, f'group_{gp}') |
| Hp_raw, Sp_raw, _, _, _ = load_structure(gp_dir, ham_file) |
| Hp_R, Sp_R = build_HR(Hp_raw, Sp_raw, site_norbits, norbits, csum) |
| eigs_p, evecs_p = compute_eigenvectors(Hp_R, Sp_R, KGRID, norbits) |
| ep_arr[d] = np.array(eigs_p) |
|
|
| delta_p = _get_delta_cart_ang(scf0_dir, gp_dir, atom) |
| S_cross_p = build_cross_overlap_R_exact(S0_raw, Sp_raw, atom, |
| site_norbits, norbits, csum, delta_p) |
| for ik in range(NK): |
| Sk_cross = build_Sk_cross(S_cross_p, KGRID[ik], norbits) |
| U_arr[d, ik, ik] = compute_brackets(evecs_p[ik], evecs_0[ik], Sk_cross) |
|
|
| |
| gm_dir = os.path.join(DISP_DIR, f'group_{gm}') |
| Hm_raw, Sm_raw, _, _, _ = load_structure(gm_dir, ham_file) |
| Hm_R, Sm_R = build_HR(Hm_raw, Sm_raw, site_norbits, norbits, csum) |
| eigs_m, evecs_m = compute_eigenvectors(Hm_R, Sm_R, KGRID, norbits) |
| epm_arr[d] = np.array(eigs_m) |
|
|
| delta_m = _get_delta_cart_ang(scf0_dir, gm_dir, atom) |
| S_cross_m = build_cross_overlap_R_exact(S0_raw, Sm_raw, atom, |
| site_norbits, norbits, csum, delta_m) |
| for ik in range(NK): |
| Sk_cross = build_Sk_cross(S_cross_m, KGRID[ik], norbits) |
| V_arr[d, ik, ik] = compute_brackets(evecs_m[ik], evecs_0[ik], Sk_cross) |
|
|
| print(f' d={d} done') |
|
|
| np.savez(out_file, |
| U_list=U_arr, V_list=V_arr, |
| ek_list=ek_arr, ep_list=ep_arr, epm_list=epm_arr) |
| print(f'\n Saved: {out_file}') |
| print(f' U_list shape: {U_arr.shape} (n_dirs, nk, nk, nbands, nbands)') |
|
|
|
|
| def run_all_brackets(): |
| print('=' * 60) |
| print('Step 3: AO bracket computation') |
| print('=' * 60) |
| for source in ('hpro', 'e3'): |
| aodir_0 = os.path.join(DISP_DIR, 'scf_0', 'reconstruction', 'aohamiltonian') |
| ham_file = 'hamiltonians.h5' if source == 'hpro' else 'hamiltonians_pred_e3.h5' |
| if not os.path.exists(os.path.join(aodir_0, ham_file)): |
| print(f' [{source.upper()}] {ham_file} not found in scf_0 — skipping') |
| continue |
| compute_ao_brackets(source) |
|
|
|
|
| |
| |
| |
|
|
| def main(): |
| parser = argparse.ArgumentParser(description='ml_epc.py: EPC AO bracket pipeline') |
| parser.add_argument('--hpro', action='store_true', help='Run only HPRO reconstruction') |
| parser.add_argument('--e3', action='store_true', help='Run only DeepH-E3 inference') |
| parser.add_argument('--brackets', action='store_true', help='Run only AO bracket computation') |
| args = parser.parse_args() |
|
|
| run_all = not (args.hpro or args.e3 or args.brackets) |
|
|
| if run_all or args.hpro: |
| run_all_hpro() |
|
|
| if run_all or args.e3: |
| run_all_e3() |
|
|
| if run_all or args.brackets: |
| run_all_brackets() |
|
|
| print('\nml_epc.py done.') |
| print('Run diamond_ml.jl to compute EPC with AO brackets and compare with DFT.') |
|
|
|
|
| if __name__ == '__main__': |
| main() |
|
|