def certify_and_datumize(dicary: Dict[str, Union[float, Decimal, np.ndarray]], *, plump: bool = False, nat: int = None) -> Dict[str, Datum]: """Convert a QCVariable: data dictionary to QCVariable: Datum dictionary, filling in details from the glossary. Parameters ---------- dicary : dict Dictionary where keys are strings (proper case – mostly UC) and values are float-like """ calcinfo = [] for pv, var in dicary.items(): if pv in qcvardefs.keys(): doi = qcvardefs[pv].get("doi", None) if plump and isinstance(var, np.ndarray) and var.ndim == 1: var = var.reshape( eval(qcvardefs[pv]["dimension"].format(nat=nat))) calcinfo.append( Datum(pv, qcvardefs[pv]["units"], var, doi=doi, glossary=qcvardefs[pv]["glossary"])) else: defined = sorted(qcvardefs.keys()) raise ValidationError("Undefined QCvar!: {}\n{}".format( pv, "\n\t".join(defined))) return {info.label: info for info in calcinfo}
def mol_to_molecule(mol: pb.Mol) -> Molecule: """Convert mol protobuf message to Molecule Note: Should not use for returning AtomicResults objects because the AtomicResult object should be a direct superset of the AtomicInput that created it (and already contains the Molecule submitted by the user) """ if mol.units == pb.Mol.UnitType.ANGSTROM: geom_angstrom = Datum("geometry", "angstrom", array(mol.xyz)) geom_bohr = geom_angstrom.to_units("bohr") elif mol.units == pb.Mol.UnitType.BOHR: geom_bohr = array(mol.xyz) else: raise ValueError( f"Unknown Unit Type: {mol.units} for molecular geometry") return Molecule( symbols=mol.atoms, geometry=geom_bohr, molecular_multiplicity=mol.multiplicity, )
def filter_nonvib(vibinfo, remove=None): """From a dictionary of vibration Datum, remove normal coordinates. Parameters ---------- vibinfo : dict of vibration Datum Results of Hessian analysis. remove : list of int, optional 0-indexed indices of normal modes to remove from `vibinfo`. If None, non-vibrations (R, T, or TR as labeled in `vibinfo['TRV']`) will be removed. Returns ------- dict of vibration Datum Copy of input `vibinfo` with the specified modes removed from all dictionary entries. Examples -------- # after a harmonic analysis, remove first translations and rotations and then all non-A1 vibs >>> allnormco = harmonic_analysis(...) >>> allvibs = filter_nonvib(allnormco) >>> a1vibs = filter_nonvib(allvibs, remove=[i for i, d in enumerate(allvibs['gamma'].data) if d != 'A1']) """ work = {} if remove is None: remove = [ idx for idx, dat in enumerate(vibinfo['TRV'].data) if dat != 'V' ] for asp, oasp in vibinfo.items(): if asp in ['q', 'w', 'x']: axis = 1 else: axis = 0 work[asp] = Datum(oasp.label, oasp.units, np.delete(oasp.data, remove, axis=axis), comment=oasp.comment, numeric=False) return work
def certify(dicary, plump=False, nat=None): """ Parameters ---------- dicary : dict Dictionary where keys are strings (proper case – mostly UC) and values are float-like """ calcinfo = [] for pv, var in dicary.items(): if pv in qcvardefs.keys(): doi = qcvardefs[pv].get('doi', None) if plump and isinstance(var, np.ndarray) and var.ndim == 1: var = var.reshape(eval(qcvardefs[pv]['dimension'].format(nat=nat))) calcinfo.append(Datum(pv, qcvardefs[pv]['units'], var, doi=doi, glossary=qcvardefs[pv]['glossary'])) else: defined = sorted(qcvardefs.keys()) raise ValidationError('Undefined QCvar!: {}\n{}'.format(pv, '\n\t'.join(defined))) return {info.label: info for info in calcinfo}
def thermo(vibinfo, T, P, multiplicity, molecular_mass, E0, sigma, rot_const, rotor_type=None): """Perform thermochemical analysis from vibrational output. Parameters ---------- E0 : float Electronic energy [Eh] at well bottom at 0 [K], :qcvar:`CURRENT ENERGY`. molecular_mass : float Mass in [u] of molecule under analysis. multiplicity: int Spin multiplicity of molecule under analysis. rot_const : ndarray of floats (3,) rotational constants in [cm^-1] of molecule under analysis. sigma : int The rotational or external symmetry number determined from the point group. rotor_type : str The rotor type for rotational stat mech purposes: RT_ATOM, RT_LINEAR, other. T : float Temperature in [K]. Psi default 298.15. Note that 273.15 is IUPAC STP. P : float Pressure in [Pa]. Psi default 101325. Note that 100000 is IUPAC STP. Returns ------- dict of Datum, str First is every thermochemistry component in atomic units along with input conditions. Second is formatted presentation of analysis. """ sm = collections.defaultdict(float) # conditions therminfo = {} therminfo['E0'] = Datum('E0', 'Eh', E0) therminfo['B'] = Datum('rotational constants', 'cm^-1', rot_const) therminfo['sigma'] = Datum('external symmetry number', '', sigma) therminfo['T'] = Datum('temperature', 'K', T) therminfo['P'] = Datum('pressure', 'Pa', P) # electronic q_elec = multiplicity sm[('S', 'elec')] = math.log(q_elec) # translational beta = 1 / (qcel.constants.kb * T) q_trans = (2.0 * np.pi * molecular_mass * qcel.constants.amu2kg / (beta * qcel.constants.h * qcel.constants.h))**1.5 * qcel.constants.na / (beta * P) sm[('S', 'trans')] = 5 / 2 + math.log(q_trans / qcel.constants.na) sm[('Cv', 'trans')] = 3 / 2 sm[('Cp', 'trans')] = 5 / 2 sm[('E', 'trans')] = 3 / 2 * T sm[('H', 'trans')] = 5 / 2 * T # rotational if rotor_type == "RT_ATOM": pass elif rotor_type == "RT_LINEAR": q_rot = 1. / (beta * sigma * 100 * qcel.constants.c * qcel.constants.h * rot_const[1]) sm[('S', 'rot')] = 1.0 + math.log(q_rot) sm[('Cv', 'rot')] = 1 sm[('Cp', 'rot')] = 1 sm[('E', 'rot')] = T else: phi_A, phi_B, phi_C = rot_const * 100 * qcel.constants.c * qcel.constants.h / qcel.constants.kb q_rot = math.sqrt( math.pi) * T**1.5 / (sigma * math.sqrt(phi_A * phi_B * phi_C)) sm[('S', 'rot')] = 3 / 2 + math.log(q_rot) sm[('Cv', 'rot')] = 3 / 2 sm[('Cp', 'rot')] = 3 / 2 sm[('E', 'rot')] = 3 / 2 * T sm[('H', 'rot')] = sm[('E', 'rot')] # vibrational vibonly = filter_nonvib(vibinfo) ZPE_cm_1 = 1 / 2 * np.sum(vibonly['omega'].data.real) omega_str = _format_omega(vibonly['omega'].data, decimals=4) imagfreqidx = np.where( vibonly['omega'].data.imag > vibonly['omega'].data.real)[0] if len(imagfreqidx): print( "Warning: thermodynamics relations excluded imaginary frequencies: {}" .format(omega_str[imagfreqidx])) filtered_theta_vib = np.delete(vibonly['theta_vib'].data, imagfreqidx, None) filtered_omega_str = np.delete(omega_str, imagfreqidx, None) rT = filtered_theta_vib / T # reduced temperature lowfreqidx = np.where(filtered_theta_vib < 900.)[0] if len(lowfreqidx): print( "Warning: used thermodynamics relations inappropriate for low-frequency modes: {}" .format(filtered_omega_str[lowfreqidx])) sm[('S', 'vib')] = np.sum(rT / np.expm1(rT) - np.log(1 - np.exp(-rT))) sm[('Cv', 'vib')] = np.sum(np.exp(rT) * (rT / np.expm1(rT))**2) sm[('Cp', 'vib')] = sm[('Cv', 'vib')] sm[('ZPE', 'vib')] = np.sum(rT) * T / 2 sm[('E', 'vib')] = sm[('ZPE', 'vib')] + np.sum(rT * T / np.expm1(rT)) sm[('H', 'vib')] = sm[('E', 'vib')] assert (abs(ZPE_cm_1 - sm[('ZPE', 'vib')] * qcel.constants.R * qcel.constants.hartree2wavenumbers * 0.001 / qcel.constants.hartree2kJmol) < 0.1) #real_vibs = np.ma.masked_where(vibinfo['omega'].data.imag > vibinfo['omega'].data.real, vibinfo['omega'].data) # compute Gibbs for term in ['elec', 'trans', 'rot', 'vib']: sm[('G', term)] = sm[('H', term)] - T * sm[('S', term)] # convert to atomic units for term in ['elec', 'trans', 'rot', 'vib']: # terms above are unitless (S, Cv, Cp) or in units of temperature (ZPE, E, H, G) as expressions are divided by R. # R [Eh/K], computed as below, slightly diff in 7th sigfig from 3.1668114e-6 (k_B in [Eh/K]) # value listed https://en.wikipedia.org/wiki/Boltzmann_constant uconv_R_EhK = qcel.constants.R / qcel.constants.hartree2kJmol for piece in ['S', 'Cv', 'Cp']: sm[(piece, term)] *= uconv_R_EhK # [mEh/K] <-- [] for piece in ['ZPE', 'E', 'H', 'G']: sm[(piece, term)] *= uconv_R_EhK * 0.001 # [Eh] <-- [K] # sum corrections and totals for piece in ['S', 'Cv', 'Cp']: for term in ['elec', 'trans', 'rot', 'vib']: sm[(piece, 'tot')] += sm[(piece, term)] for piece in ['ZPE', 'E', 'H', 'G']: for term in ['elec', 'trans', 'rot', 'vib']: sm[(piece, 'corr')] += sm[(piece, term)] sm[(piece, 'tot')] = E0 + sm[(piece, 'corr')] terms = collections.OrderedDict() terms['elec'] = ' Electronic' terms['trans'] = ' Translational' terms['rot'] = ' Rotational' terms['vib'] = ' Vibrational' terms['tot'] = 'Total' terms['corr'] = 'Correction' # package results for export for entry in sm: if entry[0] in ['S', 'Cv', 'Cp']: unit = 'mEh/K' elif entry[0] in ['ZPE', 'E', 'H', 'G']: unit = 'Eh' therminfo['_'.join(entry)] = Datum( terms[entry[1]].strip().lower() + ' ' + entry[0], unit, sm[entry]) # display format_S_Cv_Cp = """\n {:19} {:11.3f} [cal/(mol K)] {:11.3f} [J/(mol K)] {:15.8f} [mEh/K]""" format_ZPE_E_H_G = """\n {:19} {:11.3f} [kcal/mol] {:11.3f} [kJ/mol] {:15.8f} [Eh]""" uconv = np.asarray( [qcel.constants.hartree2kcalmol, qcel.constants.hartree2kJmol, 1.]) # TODO rot_const, rotor_type text = '' text += """\n ==> Thermochemistry Components <==""" text += """\n\n Entropy, S""" for term in terms: text += format_S_Cv_Cp.format(terms[term] + ' S', *sm[('S', term)] * uconv) if term == 'elec': text += """ (multiplicity = {})""".format(multiplicity) elif term == 'trans': text += """ (mol. weight = {:.4f} [u], P = {:.2f} [Pa])""".format( molecular_mass, P) elif term == 'rot': text += """ (symmetry no. = {})""".format(sigma) text += """\n\n Constant volume heat capacity, Cv""" for term in terms: text += format_S_Cv_Cp.format(terms[term] + ' Cv', *sm[('Cv', term)] * uconv) text += """\n\n Constant pressure heat capacity, Cp""" for term in terms: text += format_S_Cv_Cp.format(terms[term] + ' Cp', *sm[('Cp', term)] * uconv) del terms['tot'] terms['corr'] = 'Correction' text += """\n\n ==> Thermochemistry Energy Analysis <==""" text += """\n\n Raw electronic energy, E0""" text += """\n Total E0, Electronic energy at well bottom at 0 [K] {:15.8f} [Eh]""".format( E0) text += """\n\n Zero-point energy, ZPE_vib = Sum_i nu_i / 2""" for term in terms: text += format_ZPE_E_H_G.format(terms[term] + ' ZPE', *sm[('ZPE', term)] * uconv) if term in ['vib', 'corr']: text += """ {:15.3f} [cm^-1]""".format( sm[('ZPE', term)] * qcel.constants.hartree2wavenumbers) text += """\n Total ZPE, Electronic energy at 0 [K] {:15.8f} [Eh]""".format( sm[('ZPE', 'tot')]) text += """\n\n Thermal Energy, E (includes ZPE)""" for term in terms: text += format_ZPE_E_H_G.format(terms[term] + ' E', *sm[('E', term)] * uconv) text += """\n Total E, Electronic energy at {:7.2f} [K] {:15.8f} [Eh]""".format( T, sm[('E', 'tot')]) text += """\n\n Enthalpy, H_trans = E_trans + k_B * T""" for term in terms: text += format_ZPE_E_H_G.format(terms[term] + ' H', *sm[('H', term)] * uconv) text += """\n Total H, Enthalpy at {:7.2f} [K] {:15.8f} [Eh]""".format( T, sm[('H', 'tot')]) text += """\n\n Gibbs free energy, G = H - T * S""" for term in terms: text += format_ZPE_E_H_G.format(terms[term] + ' G', *sm[('G', term)] * uconv) text += """\n Total G, Free enthalpy at {:7.2f} [K] {:15.8f} [Eh]\n""".format( T, sm[('G', 'tot')]) return therminfo, text
def harmonic_analysis(hess, geom, mass, basisset, irrep_labels, dipder=None, project_trans=True, project_rot=True): """Like so much other Psi4 goodness, originally by @andysim Parameters ---------- hess : ndarray of float (3*nat, 3*nat) non-mass-weighted Hessian in atomic units, [Eh/a0/a0]. geom : ndarray of float (nat, 3) geometry [a0] at which Hessian computed. mass : ndarray of float (nat,) atomic masses [u]. basisset : psi4.core.BasisSet Basis set object (can be dummy, e.g., STO-3G) for SALCs. irrep_labels : list of str Irreducible representation labels. dipder : ndarray of float (3, 3 * nat) dipole derivatives in atomic units, [Eh a0/u] or [(e a0/a0)^2/u] project_trans : bool, optional Idealized translations projected out of final vibrational analysis. project_rot : bool, optional Idealized rotations projected out of final vibrational analysis. Returns ------- dict, text Returns dictionary of vibration Datum objects (fields: label units data comment). Also returns text suitable for printing. .. _`table:vibaspectinfo`: +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | key | description (label & comment) | units | data (real/imaginary modes) | +===============+============================================+===========+======================================================+ | omega | frequency | cm^-1 | nd.array(ndof) complex (real/imag) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | q | normal mode, normalized mass-weighted | a0 u^1/2 | ndarray(ndof, ndof) float | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | w | normal mode, un-mass-weighted | a0 | ndarray(ndof, ndof) float | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | x | normal mode, normalized un-mass-weighted | a0 | ndarray(ndof, ndof) float | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | degeneracy | degree of degeneracy | | ndarray(ndof) int | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | TRV | translation/rotation/vibration | | ndarray(ndof) str 'TR' or 'V' or '-' for partial | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | gamma | irreducible representation | | ndarray(ndof) str irrep or None if unclassifiable | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | mu | reduced mass | u | ndarray(ndof) float (+/+) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | k | force constant | mDyne/A | ndarray(ndof) float (+/-) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | DQ0 | RMS deviation v=0 | a0 u^1/2 | ndarray(ndof) float (+/0) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | Qtp0 | Turning point v=0 | a0 u^1/2 | ndarray(ndof) float (+/0) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | Xtp0 | Turning point v=0 | a0 | ndarray(ndof) float (+/0) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | theta_vib | char temp | K | ndarray(ndof) float (+/0) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ | IR_intensity | infrared intensity | km/mol | ndarray(ndof) float (+/+) | +---------------+--------------------------------------------+-----------+------------------------------------------------------+ Examples -------- # displacement of first atom in highest energy mode >>> vibinfo['x'].data[:, -1].reshape(nat, 3)[0] # remove translations & rotations >>> vibonly = filter_nonvib(vibinfo) """ from psi4 import core if (mass.shape[0] == geom.shape[0] == (hess.shape[0] // 3) == (hess.shape[1] // 3)) and (geom.shape[1] == 3): pass else: raise ValidationError( f"""Dimension mismatch among mass ({mass.shape}), geometry ({geom.shape}), and Hessian ({hess.shape})""" ) def mat_symm_info(a, atol=1e-14, lbl='array', stol=None): symm = np.allclose(a, a.T, atol=atol) herm = np.allclose(a, a.conj().T, atol=atol) ivrt = a.shape[0] - np.linalg.matrix_rank(a, tol=stol) return """ {:32} Symmetric? {} Hermitian? {} Lin Dep Dim? {:2}""".format( lbl + ':', symm, herm, ivrt) def vec_in_space(vec, space, tol=1.0e-4): merged = np.vstack((space, vec)) u, s, v = np.linalg.svd(merged) return (s[-1] < tol) vibinfo = {} text = [] nat = len(mass) text.append("""\n\n ==> Harmonic Vibrational Analysis <==\n""") if nat == 1: nrt_expected = 3 elif np.linalg.matrix_rank(geom) == 1: nrt_expected = 5 else: nrt_expected = 6 nmwhess = hess.copy() text.append( mat_symm_info(nmwhess, lbl='non-mass-weighted Hessian') + ' (0)') # get SALC object, possibly w/o trans & rot mints = core.MintsHelper(basisset) cdsalcs = mints.cdsalcs(0xFF, project_trans, project_rot) Uh = collections.OrderedDict() for h, lbl in enumerate(irrep_labels): tmp = np.asarray(cdsalcs.matrix_irrep(h)) if tmp.size > 0: Uh[lbl] = tmp # form projector of translations and rotations space = ('T' if project_trans else '') + ('R' if project_rot else '') TRspace = _get_TR_space(mass, geom, space=space, tol=LINEAR_A_TOL) nrt = TRspace.shape[0] text.append( f' projection of translations ({project_trans}) and rotations ({project_rot}) removed {nrt} degrees of freedom ({nrt_expected})' ) P = np.identity(3 * nat) for irt in TRspace: P -= np.outer(irt, irt) text.append(mat_symm_info(P, lbl='total projector') + f' ({nrt})') # mass-weight & solve sqrtmmm = np.repeat(np.sqrt(mass), 3) sqrtmmminv = np.divide(1.0, sqrtmmm) mwhess = np.einsum('i,ij,j->ij', sqrtmmminv, nmwhess, sqrtmmminv) text.append(mat_symm_info(mwhess, lbl='mass-weighted Hessian') + ' (0)') pre_force_constant_au = np.linalg.eigvalsh(mwhess) idx = np.argsort(pre_force_constant_au) pre_force_constant_au = pre_force_constant_au[idx] uconv_cm_1 = ( np.sqrt(qcel.constants.na * qcel.constants.hartree2J * 1.0e19) / (2 * np.pi * qcel.constants.c * qcel.constants.bohr2angstroms)) pre_frequency_cm_1 = np.lib.scimath.sqrt( pre_force_constant_au) * uconv_cm_1 pre_lowfreq = np.where(np.real(pre_frequency_cm_1) < 100.0)[0] pre_lowfreq = np.append( pre_lowfreq, np.arange(nrt_expected)) # catch at least nrt modes for lf in set(pre_lowfreq): vlf = pre_frequency_cm_1[lf] if vlf.imag > vlf.real: text.append( ' pre-proj low-frequency mode: {:9.4f}i [cm^-1]'.format( vlf.imag)) else: text.append( ' pre-proj low-frequency mode: {:9.4f} [cm^-1]'.format( vlf.real)) text.append(' pre-proj all modes:' + str(_format_omega(pre_frequency_cm_1, 4))) # project & solve mwhess_proj = np.dot(P.T, mwhess).dot(P) text.append( mat_symm_info(mwhess_proj, lbl='projected mass-weighted Hessian') + f' ({nrt})') #print('projhess = ', np.array_repr(mwhess_proj)) force_constant_au, qL = np.linalg.eigh(mwhess_proj) # expected order for vibrations is steepest downhill to steepest uphill idx = np.argsort(force_constant_au) force_constant_au = force_constant_au[idx] qL = qL[:, idx] qL = _phase_cols_to_max_element(qL) vibinfo['q'] = Datum('normal mode', 'a0 u^1/2', qL, comment='normalized mass-weighted') # frequency, LAB II.17 frequency_cm_1 = np.lib.scimath.sqrt(force_constant_au) * uconv_cm_1 vibinfo['omega'] = Datum('frequency', 'cm^-1', frequency_cm_1) # degeneracies ufreq, uinv, ucts = np.unique(np.around(frequency_cm_1, 2), return_inverse=True, return_counts=True) vibinfo['degeneracy'] = Datum('degeneracy', '', ucts[uinv]) # look among the symmetry subspaces h for one to which the normco # of vib does *not* add an extra dof to the vector space active = [] irrep_classification = [] for idx, vib in enumerate(frequency_cm_1): if vec_in_space(qL[:, idx], TRspace, 1.0e-4): active.append('TR') irrep_classification.append(None) else: active.append('V') for h in Uh.keys(): if vec_in_space(qL[:, idx], Uh[h], 1.0e-4): irrep_classification.append(h) break else: irrep_classification.append(None) # catch partial Hessians if np.linalg.norm(vib) < 1.e-3: active[-1] = '-' vibinfo['TRV'] = Datum('translation/rotation/vibration', '', active, numeric=False) vibinfo['gamma'] = Datum('irreducible representation', '', irrep_classification, numeric=False) lowfreq = np.where(np.real(frequency_cm_1) < 100.0)[0] lowfreq = np.append(lowfreq, np.arange(nrt_expected)) # catch at least nrt modes for lf in set(lowfreq): vlf = frequency_cm_1[lf] if vlf.imag > vlf.real: text.append( ' post-proj low-frequency mode: {:9.4f}i [cm^-1] ({})'.format( vlf.imag, active[lf])) else: text.append( ' post-proj low-frequency mode: {:9.4f} [cm^-1] ({})'.format( vlf.real, active[lf])) text.append(' post-proj all modes:' + str(_format_omega(frequency_cm_1, 4)) + '\n') if project_trans and not project_rot: text.append( f' Note that "Vibration"s include {nrt_expected - 3} un-projected rotation-like modes.' ) elif not project_trans and not project_rot: text.append( f' Note that "Vibration"s include {nrt_expected} un-projected rotation-like and translation-like modes.' ) # general conversion factors, LAB II.11 uconv_K = (qcel.constants.h * qcel.constants.na * 1.0e21) / (8 * np.pi * np.pi * qcel.constants.c) uconv_S = np.sqrt( (qcel.constants.c * (2 * np.pi * qcel.constants.bohr2angstroms)**2) / (qcel.constants.h * qcel.constants.na * 1.0e21)) # normco & reduced mass, LAB II.14 & II.15 wL = np.einsum('i,ij->ij', sqrtmmminv, qL) vibinfo['w'] = Datum('normal mode', 'a0', wL, comment='un-mass-weighted') reduced_mass_u = np.divide(1.0, np.linalg.norm(wL, axis=0)**2) vibinfo['mu'] = Datum('reduced mass', 'u', reduced_mass_u) xL = np.sqrt(reduced_mass_u) * wL vibinfo['x'] = Datum('normal mode', 'a0', xL, comment='normalized un-mass-weighted') # IR intensities, CCQC Proj. Eqns. 15-16 uconv_kmmol = (qcel.constants.get("Avogadro constant") * np.pi * 1.e-3 * qcel.constants.get("electron mass in u") * qcel.constants.get("fine-structure constant")**2 * qcel.constants.get("atomic unit of length") / 3) uconv_D2A2u = (qcel.constants.get('atomic unit of electric dipole mom.') * 1.e11 / qcel.constants.get('hertz-inverse meter relationship') / qcel.constants.get('atomic unit of length'))**2 if not (dipder is None or np.array(dipder).size == 0): qDD = dipder.dot(wL) ir_intensity = np.zeros(qDD.shape[1]) for i in range(qDD.shape[1]): ir_intensity[i] = qDD[:, i].dot(qDD[:, i]) # working but not needed #vibinfo['IR_intensity'] = Datum('infrared intensity', 'Eh a0/u', ir_intensity) #ir_intensity_D2A2u = ir_intensity * uconv_D2A2u #vibinfo['IR_intensity'] = Datum('infrared intensity', '(D/AA)^2/u', ir_intens_D2A2u) ir_intensity_kmmol = ir_intensity * uconv_kmmol vibinfo['IR_intensity'] = Datum('infrared intensity', 'km/mol', ir_intensity_kmmol) # force constants, LAB II.16 (real compensates for earlier sqrt) uconv_mdyne_a = (0.1 * (2 * np.pi * qcel.constants.c)**2) / qcel.constants.na force_constant_mdyne_a = reduced_mass_u * ( frequency_cm_1 * frequency_cm_1).real * uconv_mdyne_a vibinfo['k'] = Datum('force constant', 'mDyne/A', force_constant_mdyne_a) force_constant_cm_1_bb = reduced_mass_u * ( frequency_cm_1 * frequency_cm_1).real * uconv_S * uconv_S Datum('force constant', 'cm^-1/a0^2', force_constant_cm_1_bb, comment="Hooke's Law") # turning points, LAB II.20 (real & zero since turning point silly for imag modes) nu = 0 turning_point_rnc = np.sqrt(2.0 * nu + 1.0) with np.errstate(divide='ignore'): turning_point_bohr_u = turning_point_rnc / ( np.sqrt(frequency_cm_1.real) * uconv_S) turning_point_bohr_u[turning_point_bohr_u == np.inf] = 0. vibinfo['Qtp0'] = Datum('Turning point v=0', 'a0 u^1/2', turning_point_bohr_u) with np.errstate(divide='ignore'): turning_point_bohr = turning_point_rnc / ( np.sqrt(frequency_cm_1.real * reduced_mass_u) * uconv_S) turning_point_bohr[turning_point_bohr == np.inf] = 0. vibinfo['Xtp0'] = Datum('Turning point v=0', 'a0', turning_point_bohr) rms_deviation_bohr_u = turning_point_bohr_u / np.sqrt(2.0) vibinfo['DQ0'] = Datum('RMS deviation v=0', 'a0 u^1/2', rms_deviation_bohr_u) # characteristic vibrational temperature, RAK thermo & https://en.wikipedia.org/wiki/Vibrational_temperature # (imag freq zeroed) uconv_K = 100 * qcel.constants.h * qcel.constants.c / qcel.constants.kb vib_temperature_K = frequency_cm_1.real * uconv_K vibinfo['theta_vib'] = Datum('char temp', 'K', vib_temperature_K) return vibinfo, '\n'.join(text)
def nbody_gufunc(func: Union[str, Callable], method_string: str, **kwargs): """ Computes the nbody interaction energy, gradient, or Hessian depending on input. This is a generalized univeral function for computing interaction and total quantities. :returns: *return type of func* |w--w| The data. :returns: (*float*, :py:class:`~psi4.core.Wavefunction`) |w--w| data and wavefunction with energy/gradient/hessian set appropriately when **return_wfn** specified. :type func: Callable :param func: ``energy`` || etc. Python function that accepts method_string and a molecule. Returns a energy, gradient, or Hessian as requested. :type method_string: str :param method_string: ``'scf'`` || ``'mp2'`` || ``'ci5'`` || etc. First argument, lowercase and usually unlabeled. Indicates the computational method to be passed to func. :type molecule: :ref:`molecule <op_py_molecule>` :param molecule: ``h2o`` || etc. The target molecule, if not the last molecule defined. :type return_wfn: :ref:`boolean <op_py_boolean>` :param return_wfn: ``'on'`` || |dl| ``'off'`` |dr| Indicate to additionally return the :py:class:`~psi4.core.Wavefunction` calculation result as the second element of a tuple. :type bsse_type: str or list :param bsse_type: ``'cp'`` || ``['nocp', 'vmfc']`` || |dl| ``None`` |dr| || etc. Type of BSSE correction to compute: CP, NoCP, or VMFC. The first in this list is returned by this function. By default, this function is not called. :type max_nbody: int :param max_nbody: ``3`` || etc. Maximum n-body to compute, cannot exceed the number of fragments in the moleucle. :type ptype: str :param ptype: ``'energy'`` || ``'gradient'`` || ``'hessian'`` Type of the procedure passed in. :type return_total_data: :ref:`boolean <op_py_boolean>` :param return_total_data: ``'on'`` || |dl| ``'off'`` |dr| If True returns the total data (energy/gradient/etc) of the system, otherwise returns interaction data. :type levels: dict :param levels: ``{1: 'ccsd(t)', 2: 'mp2', 'supersystem': 'scf'}`` || ``{1: 2, 2: 'ccsd(t)', 3: 'mp2'}`` || etc Dictionary of different levels of theory for different levels of expansion Note that method_string is not used in this case. supersystem computes all higher order n-body effects up to nfragments. :type embedding_charges: dict :param embedding_charges: ``{1: [-0.834, 0.417, 0.417], ..}`` Dictionary of atom-centered point charges. keys: 1-based index of fragment, values: list of charges for each fragment. :type charge_method: str :param charge_method: ``scf/6-31g`` || ``b3lyp/6-31g*`` || etc Method to compute point charges for monomers. Overridden by embedding_charges if both are provided. :type charge_type: str :param charge_type: ``MULLIKEN_CHARGES`` || ``LOWDIN_CHARGES`` Default is ``MULLIKEN_CHARGES`` """ # Initialize dictionaries for easy data passing metadata, component_results, nbody_results = {}, {}, {} # Parse some kwargs kwargs = driver_util.kwargs_lower(kwargs) if kwargs.get("levels", False): return driver_nbody_helper.multi_level(func, **kwargs) metadata["ptype"] = kwargs.pop("ptype", None) metadata["return_wfn"] = kwargs.pop("return_wfn", False) metadata["return_total_data"] = kwargs.pop("return_total_data", False) metadata["molecule"] = kwargs.pop("molecule") metadata["molecule"].update_geometry() metadata["molecule"].fix_com(True) metadata["molecule"].fix_orientation(True) metadata["molecule"].update_geometry() # MM metadata["embedding_charges"] = kwargs.get("embedding_charges", False) metadata["kwargs"] = kwargs # core.clean_variables() if metadata["ptype"] not in ["energy", "gradient", "hessian"]: raise ValidationError( """N-Body driver: The ptype '%s' is not recognized.""" % metadata["ptype"]) if metadata["ptype"] not in ["energy"]: raise ValidationError( """N-Body driver: Only energy ptype available for now.""") # Parse bsse_type, raise exception if not provided or unrecognized metadata["bsse_type_list"] = kwargs.pop("bsse_type") if metadata["bsse_type_list"] is None: raise ValidationError("N-Body GUFunc: Must pass a bsse_type") if not isinstance(metadata["bsse_type_list"], list): metadata["bsse_type_list"] = [metadata["bsse_type_list"]] for num, btype in enumerate(metadata["bsse_type_list"]): metadata["bsse_type_list"][num] = btype.lower() if btype.lower() not in ["cp", "nocp", "vmfc"]: raise ValidationError( "N-Body GUFunc: bsse_type '%s' is not recognized" % btype.lower()) metadata["max_nbody"] = kwargs.get("max_nbody", -1) if metadata["molecule"].nfragments() == 1: raise ValidationError( "N-Body requires active molecule to have more than 1 fragment.") metadata["max_frag"] = metadata["molecule"].nfragments() if metadata["max_nbody"] == -1: metadata["max_nbody"] = metadata["molecule"].nfragments() else: metadata["max_nbody"] = min(metadata["max_nbody"], metadata["max_frag"]) # Flip this off for now, needs more testing # If we are doing CP lets save them integrals # if 'cp' in bsse_type_list and (len(bsse_type_list) == 1): # # Set to save RI integrals for repeated full-basis computations # ri_ints_io = core.get_global_option('DF_INTS_IO') # # inquire if above at all applies to dfmp2 or just scf # core.set_global_option('DF_INTS_IO', 'SAVE') # psioh = core.IOManager.shared_object() # psioh.set_specific_retention(97, True) bsse_str = metadata["bsse_type_list"][0] if len(metadata["bsse_type_list"]) > 1: bsse_str = str(metadata["bsse_type_list"]) print("\n\n") print(" ===> N-Body Interaction Abacus <===\n") print(" BSSE Treatment: %s\n" % bsse_str) # Get compute list metadata = build_nbody_compute_list(metadata) # Compute N-Body components component_results = compute_nbody_components(func, method_string, metadata) # Assemble N-Body quantities nbody_results = assemble_nbody_components(metadata, component_results) # Build wfn and bind variables jobrec = {} calcinfo = [] # TODO hack jobrec["molecule"] = metadata["molecule"] # wfn = core.Wavefunction.build(metadata['molecule'], 'def2-svp') dicts = [ "energies", "ptype", "intermediates", "energy_body_dict", "gradient_body_dict", "hessian_body_dict", "nbody", "cp_energy_body_dict", "nocp_energy_body_dict", "vmfc_energy_body_dict", ] if metadata["ptype"] == "gradient": # wfn.set_gradient(nbody_results['ret_ptype']) calcinfo.append( Datum("CURRENT GRADIENT", "Eh/a0", nbody_results["ret_ptype"])) nbody_results["gradient_body_dict"] = nbody_results["ptype_body_dict"] elif metadata["ptype"] == "hessian": nbody_results["hessian_body_dict"] = nbody_results["ptype_body_dict"] # wfn.set_hessian(nbody_results['ret_ptype']) calcinfo.append( Datum("CURRENT HESSIAN", "Eh/a0/a0", nbody_results["ret_ptype"])) component_results_gradient = component_results.copy() component_results_gradient["ptype"] = component_results_gradient[ "gradients"] metadata["ptype"] = "gradient" nbody_results_gradient = assemble_nbody_components( metadata, component_results_gradient) # wfn.set_gradient(nbody_results_gradient['ret_ptype']) nbody_results["gradient_body_dict"] = nbody_results_gradient[ "ptype_body_dict"] calcinfo.append( Datum("CURRENT GRADIENT", "Eh/a0", nbody_results_gradient["ret_ptype"])) for r in [component_results, nbody_results]: for d in r: if d in dicts: for var, value in r[d].items(): try: # wfn.set_scalar_variable(str(var), value) # core.set_scalar_variable(str(var), value) calcinfo.append(Datum(str(var), "Eh", value)) except: # wfn.set_array_variable(d.split('_')[0].upper() + ' ' + str(var), core.Matrix.from_array(value)) calcinfo.append( Datum( d.split("_")[0].upper() + " " + str(var), "?", value)) # core.set_variable("CURRENT ENERGY", nbody_results['ret_energy']) # wfn.set_variable("CURRENT ENERGY", nbody_results['ret_energy']) calcinfo.append(Datum("CURRENT ENERGY", "Eh", nbody_results["ret_energy"])) if metadata["ptype"] == "gradient": # core.set_variable("CURRENT GRADIENT", nbody_results['ret_ptype']) calcinfo.append( Datum("CURRENT GRADIENT", "Eh/a0", nbody_results["ret_ptype"])) elif metadata["ptype"] == "hessian": # core.set_variable("CURRENT HESSIAN", nbody_results['ret_ptype']) calcinfo.append( Datum("CURRENT HESSIAN", "Eh/a0/a0", nbody_results["ret_ptype"])) jobrec["qcvars"] = {info.label: info for info in calcinfo} pp.pprint(jobrec) pe.active_qcvars = copy.deepcopy(jobrec["qcvars"]) if metadata["return_wfn"]: return (nbody_results["ret_ptype"], jobrec) else: return nbody_results["ret_ptype"]