octadist.io

octadist.src.io.is_cif(f)[source]

Check if the input file is .cif file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is CIF file, return True.

See also

get_coord_cif

Find atomic coordinates of molecule from CIF file.

Notes

More details about CIF file format are provided at https://en.wikipedia.org/wiki/Crystallographic_Information_File.

Examples

>>> # example.cif
>>> # example
>>> # _audit_creation_date              2012-10-26T21:09:50-0400
>>> # _audit_creation_method            fapswitch 2.2
>>> # _symmetry_space_group_name_H-M    P1
>>> # _symmetry_Int_Tables_number       1
>>> # _space_group_crystal_system       triclinic
>>> # _cell_length_a                    16.012374
>>> # _cell_length_b                    14.740457
>>> # _cell_length_c                    19.436146
>>> # _cell_angle_alpha                 89.939227
>>> # _cell_angle_beta                  90.110039
>>> # _cell_angle_gamma                 90.015104
>>> # _cell_volume                      4587.49671393
>>> #
>>> # loop_
>>> # _atom_site_label
>>> # _atom_site_type_symbol
>>> # _atom_type_description
>>> # _atom_site_fract_x
>>> # _atom_site_fract_y
>>> # _atom_site_fract_z
>>> # _atom_type_partial_charge
>>> # C1    C     C_R   0.340882 0.499989 0.500098 0.541130
>>> # C2    C     C_R   0.528123 0.048033 0.558069 0.232589
>>> # C3    C     C_R   0.499931 0.902862 0.500001 -0.063750
>>> # C4    C     C_R   0.500061 0.097137 0.500001 -0.063745
>>> # C5    C     C_1   0.499958 0.802655 0.499991 0.266033
>>> # ...
>>> is_cif("example.cif")
True
octadist.src.io.get_coord_cif(f)[source]

Get coordinate from .cif file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "example.cif"
>>> atom, coord = get_coord_cif(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.is_xyz(f)[source]

Check if the input file is .xyz file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is XYZ file, return True.

See also

get_coord_xyz

Find atomic coordinates of molecule from XYZ file.

Examples

>>> # example.xyz
>>> # 20
>>> # Comment: From Excel file
>>> # Fe  6.251705    9.063211    5.914842
>>> # N   8.15961     9.066456    5.463087
>>> # N   6.749414    10.457551   7.179682
>>> # N   5.709997    10.492955   4.658257
>>> # N   4.350474    9.106286    6.356091
>>> # O   5.789096    7.796326    4.611355
>>> # O   6.686381    7.763872    7.209699
>>> # ...
>>> is_xyz("example.xyz")
True
octadist.src.io.get_coord_xyz(f)[source]

Get coordinate from .xyz file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "Fe-distorted-complex.xyz"
>>> atom, coord = get_coord_xyz(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.is_gaussian(f)[source]

Check if the input file is Gaussian file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is Gaussian output file, return True.

See also

get_coord_gaussian

Find atomic coordinates of molecule from Gaussian file.

Examples

>>> # gaussian.log
>>> #                             Standard orientation:
>>> # ---------------------------------------------------------------------
>>> # Center     Atomic      Atomic             Coordinates (Angstroms)
>>> # Number     Number       Type             X           Y           Z
>>> # ---------------------------------------------------------------------
>>> #      1         26           0        0.000163    1.364285   -0.000039
>>> #      2          8           0        0.684192    0.084335   -1.192008
>>> #      3          8           0       -0.683180    0.083251    1.191173
>>> #      4          7           0        1.639959    1.353157    1.006941
>>> #      5          7           0       -0.563377    2.891083    1.435925
>>> # ...
>>> is_gaussian("gaussian.log")
True
octadist.src.io.get_coord_gaussian(f)[source]

Extract XYZ coordinate from Gaussian output file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "Gaussian-Fe-distorted-complex.out"
>>> atom, coord = get_coord_gaussian(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.is_nwchem(f)[source]

Check if the input file is NWChem file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is NWChem output file, return True.

See also

get_coord_nwchem

Find atomic coordinates of molecule from NWChem file.

Examples

>>> # nwchem.out
>>> #   ----------------------
>>> #   Optimization converged
>>> #   ----------------------
>>> # ...
>>> # ...
>>> #  No.       Tag          Charge          X              Y              Z
>>> # ---- ---------------- ---------- -------------- -------------- --------------
>>> #    1 Ru(Fragment=1)      44.0000    -3.04059115    -0.08558108    -0.07699482
>>> #    2 C(Fragment=1)        6.0000    -1.62704660     2.40971357     0.63980357
>>> #    3 C(Fragment=1)        6.0000    -0.61467778     0.59634595     1.68841986
>>> #    4 C(Fragment=1)        6.0000     0.31519183     1.41684566     2.30745116
>>> #    5 C(Fragment=1)        6.0000     0.28773462     2.80126911     2.08006241
>>> # ...
>>> is_nwchem("nwchem.out")
True
octadist.src.io.get_coord_nwchem(f)[source]

Extract XYZ coordinate from NWChem output file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "NWChem-Fe-distorted-complex.out"
>>> atom, coord = get_coord_nwchem(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.is_orca(f)[source]

Check if the input file is ORCA file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is ORCA output file, return True.

See also

get_coord_orca

Find atomic coordinates of molecule from ORCA file.

Examples

>>> # orca.out
>>> # ---------------------------------
>>> # CARTESIAN COORDINATES (ANGSTROEM)
>>> # ---------------------------------
>>> #   C      0.009657    0.000000    0.005576
>>> #   C      0.009657   -0.000000    1.394424
>>> #   C      1.212436   -0.000000    2.088849
>>> #   C      2.415214    0.000000    1.394425
>>> #   C      2.415214   -0.000000    0.005575
>>> # ...
>>> is_orca("orca.out")
True
octadist.src.io.get_coord_orca(f)[source]

Extract XYZ coordinate from ORCA output file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "ORCA-Fe-distorted-complex.out"
>>> atom, coord = get_coord_orca(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.is_qchem(f)[source]

Check if the input file is Q-Chem file format.

Parameters:

f (str) – User input filename.

Returns:

bool (bool) – If file is Q-Chem output file, return True.

See also

get_coord_qchem

Find atomic coordinates of molecule from Q-Chem file.

Examples

>>> # qchem.out
>>> # ******************************
>>> # **  OPTIMIZATION CONVERGED  **
>>> # ******************************
>>> #                            Coordinates (Angstroms)
>>> #     ATOM                X               Y               Z
>>> #      1  C         0.2681746845   -0.8206222796   -0.3704019386
>>> #      2  C        -1.1809302341   -0.5901746612   -0.6772716414
>>> #      3  H        -1.6636318262   -1.5373167851   -0.9496501352
>>> #      4  H        -1.2829834971    0.0829227646   -1.5389938241
>>> #      5  C        -1.9678565203    0.0191922768    0.5346693165
>>> # ...
>>> is_qchem("qchem.out")
True
octadist.src.io.get_coord_qchem(f)[source]

Extract XYZ coordinate from Q-Chem output file.

Parameters:

f (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

Examples

>>> file = "Qchem-Fe-distorted-complex.out"
>>> atom, coord = get_coord_qchem(file)
>>> atom
['Fe', 'O', 'O', 'N', 'N', 'N', 'N']
>>> coord
array([[18.268051, 11.28912 ,  2.565804],
       [19.823874, 10.436314,  1.381569],
       [19.074466,  9.706294,  3.743576],
       [17.364238, 10.733354,  0.657318],
       [16.149538, 11.306661,  2.913619],
       [18.599941, 12.116308,  4.528988],
       [18.364987, 13.407634,  2.249608]])
octadist.src.io.count_line(file=None)[source]

Count lines in an input file.

Parameters:

file (str) – Absolute or full path of input file.

Returns:

i + 1 (int) – Number of line in file.

Examples

>>> file = "[Fe(1-bpp)2][BF4]2-HS.xyz"
>>> count_line(file)
27
octadist.src.io.extract_coord(file=None)[source]

Check file type, read data, extract atomic symbols and cartesian coordinate from a structure input file provided by the user. This function can efficiently manupulate I/O process. File types currently supported are listed in notes below. Other file formats can also be implemented easily within this module.

Parameters:

file (str) – User input filename.

Returns:

  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

See also

octadist.main.OctaDist.open_file

Open file dialog and to browse input file.

octadist.main.OctaDist.search_coord

Search octahedral structure in complex.

Notes

The following are file types supported by the current virsion of OctaDist:

  • CIF

  • XYZ

  • Gaussian

  • NWChem

  • ORCA

  • Q-Chem

Examples

>>> file = "[Fe(1-bpp)2][BF4]2-HS.xyz"
>>> atom, coord = extract_coord(file)
>>> atom
['Fe', 'N', 'N', 'N', 'N', 'N', 'N', 'C', 'C']
>>> coord
array([[-1.95348286e+00,  4.51770478e+00,  1.47855811e+01],
       [-1.87618286e+00,  4.48070478e+00,  1.26484811e+01],
       [-2.18698286e+00,  4.34540478e+00,  1.69060811e+01],
       [-4.88286000e-03,  3.69060478e+00,  1.42392811e+01],
       [-1.17538286e+00,  6.38340478e+00,  1.56457811e+01],
       [-2.75078286e+00,  2.50260478e+00,  1.51806811e+01],
       [-3.90128286e+00,  5.27750478e+00,  1.40814811e+01],
       [-6.14953418e+00,  8.30666180e+00,  2.91361978e+01],
       [-8.59995241e+00,  7.11630815e+00,  4.52898814e+01]])
octadist.src.io.find_metal(atom=None, coord=None)[source]

Count the number of metal center atom in complex.

Parameters:
  • atom (list or None) – Full atomic labels of complex. Default is None.

  • coord (array_like or None) – Full atomic coordinates of complex. Default is None.

Returns:

  • atom_metal (list) – Atomic labels of metal center atom.

  • coord_metal (array_like) – Atomic coordinates of metal center atom.

  • index_metal (list) – Indices of metal atoms found.

See also

octadist.src.elements.check_atom

Convert atomic number to atomic symbol and vice versa.

Examples

>>> atom = ['Fe', 'N', 'N', 'N', 'N', 'N', 'N']
>>> coord = [[-1.95348286e+00,  4.51770478e+00,  1.47855811e+01],
             [-1.87618286e+00,  4.48070478e+00,  1.26484811e+01],
             [-3.90128286e+00,  5.27750478e+00,  1.40814811e+01],
             [-4.88286000e-03,  3.69060478e+00,  1.42392811e+01],
             [-2.18698286e+00,  4.34540478e+00,  1.69060811e+01],
             [-1.17538286e+00,  6.38340478e+00,  1.56457811e+01],
             [-2.75078286e+00,  2.50260478e+00,  1.51806811e+01]]
>>> atom_metal, coord_metal = find_metal(atom, coord)
>>> atom_metal
['Fe']
>>> coord_metal
array([[-1.95348286,  4.51770478, 14.7855811 ]])
octadist.src.io.extract_octa(atom, coord, ref_index=0, cutoff_ref_ligand=2.8)[source]

Search the octahedral structure in complex and return atoms and coordinates.

Parameters:
  • atom (list) – Full atomic labels of complex.

  • coord (array_like) – Full atomic coordinates of complex.

  • ref_index (int) – Index of the reference to be used as the center atom for neighbor atoms in octahedral structure of the complex. Python-based index. Default is 0.

  • cutoff_ref_ligand (float, optional) – Cutoff distance for screening bond distance between reference and ligand atoms. Default is 2.8.

Returns:

  • atom_octa (list) – Atomic labels of octahedral structure.

  • coord_octa (array_like) – Atomic coordinates of octahedral structure.

See also

find_metal

Find metals in complex.

octadist.main.OctaDist.search_coord

Search octahedral structure in complex.

Examples

>>> atom = ['Fe', 'N', 'N', 'N', 'N', 'N', 'N', 'C', 'C']
>>> coord = [[-1.95348286e+00,  4.51770478e+00,  1.47855811e+01],
             [-1.87618286e+00,  4.48070478e+00,  1.26484811e+01],
             [-3.90128286e+00,  5.27750478e+00,  1.40814811e+01],
             [-4.88286000e-03,  3.69060478e+00,  1.42392811e+01],
             [-2.18698286e+00,  4.34540478e+00,  1.69060811e+01],
             [-1.17538286e+00,  6.38340478e+00,  1.56457811e+01],
             [-2.75078286e+00,  2.50260478e+00,  1.51806811e+01],
             [-6.14953418e+00,  8.30666180e+00,  2.91361978e+01],
             [-8.59995241e+00,  7.11630815e+00,  4.52898814e+01]]
>>> atom_octa, coord_octa = extract_octa(atom, coord)
>>> atom_octa
['Fe', 'N', 'N', 'N', 'N', 'N', 'N']
>>> coord_octa
array([[-1.95348286e+00,  4.51770478e+00,  1.47855811e+01],
       [-1.87618286e+00,  4.48070478e+00,  1.26484811e+01],
       [-2.18698286e+00,  4.34540478e+00,  1.69060811e+01],
       [-4.88286000e-03,  3.69060478e+00,  1.42392811e+01],
       [-1.17538286e+00,  6.38340478e+00,  1.56457811e+01],
       [-2.75078286e+00,  2.50260478e+00,  1.51806811e+01],
       [-3.90128286e+00,  5.27750478e+00,  1.40814811e+01]])