Skip to content

Commit

Permalink
Merge pull request #1299 from espottesmith/master
Browse files Browse the repository at this point in the history
Removes build_MoleculeGraph, replaces with MoleculeGraph constructors
  • Loading branch information
shyuep committed Oct 9, 2018
2 parents 8f2fbdb + dd8fccd commit ee27a3b
Show file tree
Hide file tree
Showing 7 changed files with 184 additions and 131 deletions.
18 changes: 9 additions & 9 deletions pymatgen/analysis/bond_dissociation.py
Expand Up @@ -9,7 +9,7 @@
from monty.json import MSONable

from pymatgen.core.structure import Molecule
from pymatgen.analysis.graphs import build_MoleculeGraph, MolGraphSplitError
from pymatgen.analysis.graphs import MoleculeGraph, MolGraphSplitError
from pymatgen.analysis.local_env import OpenBabelNN
from pymatgen.analysis.fragmenter import open_ring
from pymatgen.io.babel import BabelMolAdaptor
Expand Down Expand Up @@ -88,10 +88,10 @@ def __init__(self, molecule_entry, fragment_entries, allow_additional_charge_sep
self.expected_charges = [molecule_entry["final_molecule"]["charge"]-2, molecule_entry["final_molecule"]["charge"]-1, molecule_entry["final_molecule"]["charge"], molecule_entry["final_molecule"]["charge"]+1]

# Build principle molecule graph
self.mol_graph = build_MoleculeGraph(Molecule.from_dict(molecule_entry["final_molecule"]),
strategy=OpenBabelNN,
reorder=False,
extend_structure=False)
self.mol_graph = MoleculeGraph.with_local_env_strategy(Molecule.from_dict(molecule_entry["final_molecule"]),
OpenBabelNN(),
reorder=False,
extend_structure=False)
# Loop through bonds, aka graph edges, and fragment and process:
for bond in self.mol_graph.graph.edges:
bonds = [(bond[0],bond[1])]
Expand Down Expand Up @@ -246,12 +246,12 @@ def filter_fragment_entries(self,fragment_entries):
elif entry["pcm_dielectric"] != self.molecule_entry["pcm_dielectric"]:
raise RuntimeError("Principle molecule has a PCM dielectric of " + str(self.molecule_entry["pcm_dielectric"]) + " but a fragment entry has a different PCM dielectric! Please only pass fragment entries with PCM details consistent with the principle entry. Exiting...")
# Build initial and final molgraphs:
entry["initial_molgraph"] = build_MoleculeGraph(Molecule.from_dict(entry["initial_molecule"]),
strategy=OpenBabelNN,
entry["initial_molgraph"] = MoleculeGraph.with_local_env_strategy(Molecule.from_dict(entry["initial_molecule"]),
OpenBabelNN(),
reorder=False,
extend_structure=False)
entry["final_molgraph"] = build_MoleculeGraph(Molecule.from_dict(entry["final_molecule"]),
strategy=OpenBabelNN,
entry["final_molgraph"] = MoleculeGraph.with_local_env_strategy(Molecule.from_dict(entry["final_molecule"]),
OpenBabelNN(),
reorder=False,
extend_structure=False)
# Classify any potential structural change that occured during optimization:
Expand Down
13 changes: 7 additions & 6 deletions pymatgen/analysis/fragmenter.py
Expand Up @@ -7,7 +7,7 @@
import logging

from monty.json import MSONable
from pymatgen.analysis.graphs import build_MoleculeGraph, MolGraphSplitError
from pymatgen.analysis.graphs import MoleculeGraph, MolGraphSplitError
from pymatgen.analysis.local_env import OpenBabelNN
from pymatgen.io.babel import BabelMolAdaptor

Expand Down Expand Up @@ -52,11 +52,12 @@ def __init__(self, molecule, edges=None, depth=1, open_rings=True, opt_steps=100
self.opt_steps = opt_steps

if edges is None:
self.mol_graph = build_MoleculeGraph(molecule, strategy=OpenBabelNN,
reorder=False, extend_structure=False)
self.mol_graph = MoleculeGraph.with_local_env_strategy(molecule, OpenBabelNN(),
reorder=False,
extend_structure=False)
else:
edges = [(e[0], e[1], {}) for e in edges]
self.mol_graph = build_MoleculeGraph(molecule, edges=edges)
edges = {(e[0], e[1]): None for e in edges}
self.mol_graph = MoleculeGraph.with_edges(molecule, edges)

self.unique_fragments = []
self.unique_fragments_from_ring_openings = []
Expand Down Expand Up @@ -161,4 +162,4 @@ def open_ring(mol_graph, bond, opt_steps):
obmol = BabelMolAdaptor.from_molecule_graph(mol_graph)
obmol.remove_bond(bond[0][0]+1, bond[0][1]+1)
obmol.localopt(steps=opt_steps)
return build_MoleculeGraph(obmol.pymatgen_mol, strategy=OpenBabelNN, reorder=False, extend_structure=False)
return MoleculeGraph.with_local_env_strategy(obmol.pymatgen_mol, OpenBabelNN(), reorder=False, extend_structure=False)
9 changes: 5 additions & 4 deletions pymatgen/analysis/functional_groups.py
Expand Up @@ -2,20 +2,20 @@
# Copyright (c) Pymatgen Development Team.
# Distributed under the terms of the MIT License.

import copy

from pymatgen.core.structure import Molecule
from pymatgen.io.babel import BabelMolAdaptor
from pymatgen.analysis.graphs import MoleculeGraph
from pymatgen.analysis.local_env import OpenBabelNN


try:
import networkx as nx
import networkx.algorithms.isomorphism as iso
except ImportError:
raise ImportError("pymatgen.analysis.functional_groups requires the "
"NetworkX graph library to be installed.")


__author__ = "Evan Spotte-Smith"
__version__ = "0.1"
__maintainer__ = "Evan Spotte-Smith"
Expand Down Expand Up @@ -291,14 +291,15 @@ def get_basic_functional_groups(self, func_groups=None):

if num_deviants <= 1:
for node in ring:
ring_group = copy.deepcopy(ring)
neighbors = self.molgraph.graph[node]

# Add hydrogens to the functional group
for neighbor in neighbors.keys():
if neighbor in hydrogens:
ring.add(neighbor)
ring_group.add(neighbor)

results.append(ring)
results.append(ring_group)

return results

Expand Down
185 changes: 116 additions & 69 deletions pymatgen/analysis/graphs.py
Expand Up @@ -136,6 +136,61 @@ def with_empty_graph(cls, structure, name="bonds",

return cls(structure, graph_data=graph_data)

@staticmethod
def with_edges(structure, edges):
"""
Constructor for MoleculeGraph, using pre-existing or pre-defined edges
with optional edge parameters.
:param molecule: Molecule object
:param edges: dict representing the bonds of the functional
group (format: {(from_index, to_index, from_image, to_image): props},
where props is a dictionary of properties, including weight.
Props should be None if no additional properties are to be
specified.
:return: sg, a StructureGraph
"""

sg = StructureGraph.with_empty_graph(structure, name="bonds",
edge_weight_name="weight",
edge_weight_units="")

for edge, props in edges.items():

try:
from_index = edge[0]
to_index = edge[1]
from_image = edge[2]
to_image = edge[3]
except TypeError:
raise ValueError("Edges must be given as (from_index, to_index,"
" from_image, to_image) tuples")

if props is not None:
if "weight" in props.keys():
weight = props["weight"]
del props["weight"]
else:
weight = None

if len(props.items()) == 0:
props = None
else:
weight = None

nodes = sg.graph.nodes
if not (from_index in nodes and to_index in nodes):
raise ValueError("Edges cannot be added if nodes are not"
" present in the graph. Please check your"
" indices.")

sg.add_edge(from_index, to_index, from_jimage=from_image,
to_jimage=to_image, weight=weight,
edge_properties=props)

sg.set_node_attributes()
return sg

@staticmethod
def with_local_env_strategy(structure, strategy):
"""
Expand Down Expand Up @@ -1454,6 +1509,8 @@ def __init__(self, molecule, graph_data=None):
if 'from_jimage' in d:
d['from_jimage'] = tuple(d['from_jimage'])

self.set_node_attributes()

@classmethod
def with_empty_graph(cls, molecule, name="bonds",
edge_weight_name=None,
Expand Down Expand Up @@ -1490,6 +1547,57 @@ def with_empty_graph(cls, molecule, name="bonds",

return cls(molecule, graph_data=graph_data)

@staticmethod
def with_edges(molecule, edges):
"""
Constructor for MoleculeGraph, using pre-existing or pre-defined edges
with optional edge parameters.
:param molecule: Molecule object
:param edges: dict representing the bonds of the functional
group (format: {(u, v): props}, where props is a dictionary of
properties, including weight. Props should be None if no
additional properties are to be specified.
:return: mg, a MoleculeGraph
"""

mg = MoleculeGraph.with_empty_graph(molecule, name="bonds",
edge_weight_name="weight",
edge_weight_units="")

for edge, props in edges.items():

try:
from_index = edge[0]
to_index = edge[1]
except TypeError:
raise ValueError("Edges must be given as (from_index, to_index)"
"tuples")

if props is not None:
if "weight" in props.keys():
weight = props["weight"]
del props["weight"]
else:
weight = None

if len(props.items()) == 0:
props = None
else:
weight = None

nodes = mg.graph.nodes
if not (from_index in nodes and to_index in nodes):
raise ValueError("Edges cannot be added if nodes are not"
" present in the graph. Please check your"
" indices.")

mg.add_edge(from_index, to_index, weight=weight,
edge_properties=props)

mg.set_node_attributes()
return mg

@staticmethod
def with_local_env_strategy(molecule, strategy, reorder=True,
extend_structure=True):
Expand Down Expand Up @@ -1557,6 +1665,7 @@ def with_local_env_strategy(molecule, strategy, reorder=True,
for duplicate in duplicates:
mg.graph.remove_edge(duplicate[0], duplicate[1], key=duplicate[2])

mg.set_node_attributes()
return mg

@property
Expand Down Expand Up @@ -1936,19 +2045,17 @@ def build_unique_fragments(self):
species = nx.get_node_attributes(remapped, "specie")
coords = nx.get_node_attributes(remapped, "coords")

edges = []
edges = {}

for from_index, to_index, key in remapped.edges:
edge_props = fragment.get_edge_data(from_index, to_index, key=key)

edges.append((from_index, to_index, edge_props))
edges[(from_index, to_index)] = edge_props

unique_mol_graphs.append(build_MoleculeGraph(Molecule(species=species,
coords=coords,
charge=self.molecule.charge),
edges=edges,
reorder=False,
extend_structure=False))
unique_mol_graphs.append(self.with_edges(Molecule(species=species,
coords=coords,
charge=self.molecule.charge),
edges))
return unique_mol_graphs

def substitute_group(self, index, func_grp, strategy, bond_order=1,
Expand Down Expand Up @@ -2666,64 +2773,4 @@ def diff(self, other, strict=True):
'other': edges_other - edges,
'both': edges.intersection(edges_other),
'dist': jaccard_dist
}


def build_MoleculeGraph(molecule, edges=None, strategy=None,
strategy_params=None, reorder=True,
extend_structure=True):
"""
General out-of-class constructor for MoleculeGraph.
:param molecule: pymatgen.core.Molecule object.
:param edges: List of tuples (from, to, properties) representing edges in
the graph. Default None.
:param strat: an instance of a
:Class: `pymatgen.analysis.local_env.NearNeighbors` object. Default
None.
:param strategy_params: dict of parameters to be passed to NearNeighbors
strategy.
:param reorder: bool, representing if graph nodes need to be reordered
following the application of the local_env strategy
:param extend_structure: If True (default), then a large artificial box
will be placed around the Molecule, because some strategies assume
periodic boundary conditions.
:return: MoleculeGraph object.
"""

if edges is None:
if strategy is not None:
if strategy_params is None:
strategy_params = {}
strat = strategy(**strategy_params)
mol_graph = MoleculeGraph.with_local_env_strategy(molecule, strat,
reorder=reorder,
extend_structure=extend_structure)
else:
raise ValueError("Must supply either edge list or"
" pymatgen.analysis.local_env strategy.")
else:
mol_graph = MoleculeGraph.with_empty_graph(molecule)
for from_index, to_index, properties in edges:
if properties is not None:
if "weight" in properties.keys():
weight = properties["weight"]
del properties["weight"]
else:
weight = None

if len(properties.items()) == 0:
properties = None
else:
weight = None

nodes = mol_graph.graph.nodes
if not (from_index in nodes and to_index in nodes):
raise ValueError("Edges cannot be added if nodes are not"
" present in the graph. Please check your"
" indices.")
mol_graph.add_edge(from_index, to_index, weight=weight,
edge_properties=properties)

mol_graph.set_node_attributes()
return mol_graph
}
10 changes: 5 additions & 5 deletions pymatgen/analysis/tests/test_fragmenter.py
Expand Up @@ -6,7 +6,7 @@
import unittest

from pymatgen.core.structure import Molecule
from pymatgen.analysis.graphs import build_MoleculeGraph
from pymatgen.analysis.graphs import MoleculeGraph
from pymatgen.analysis.local_env import OpenBabelNN
from pymatgen.util.testing import PymatgenTest
from pymatgen.analysis.fragmenter import Fragmenter
Expand Down Expand Up @@ -53,8 +53,8 @@ def test_babel_PC_defaults(self):
fragmenter = Fragmenter(molecule=self.pc)
self.assertEqual(fragmenter.open_rings,True)
self.assertEqual(fragmenter.opt_steps,10000)
default_mol_graph = build_MoleculeGraph(self.pc, strategy=OpenBabelNN,
reorder=False, extend_structure=False)
default_mol_graph = MoleculeGraph.with_local_env_strategy(self.pc, OpenBabelNN(),
reorder=False, extend_structure=False)
self.assertEqual(fragmenter.mol_graph,default_mol_graph)
self.assertEqual(len(fragmenter.unique_fragments), 13)
self.assertEqual(len(fragmenter.unique_fragments_from_ring_openings), 5)
Expand All @@ -63,8 +63,8 @@ def test_edges_given_PC_not_defaults(self):
fragmenter = Fragmenter(molecule=self.pc, edges=self.pc_edges, depth=2, open_rings=False, opt_steps=0)
self.assertEqual(fragmenter.open_rings,False)
self.assertEqual(fragmenter.opt_steps,0)
edges = [(e[0], e[1], {}) for e in self.pc_edges]
default_mol_graph = build_MoleculeGraph(self.pc, edges=edges)
edges = {(e[0], e[1]): None for e in self.pc_edges}
default_mol_graph = MoleculeGraph.with_edges(self.pc, edges=edges)
self.assertEqual(fragmenter.mol_graph,default_mol_graph)
self.assertEqual(len(fragmenter.unique_fragments), 20)
self.assertEqual(len(fragmenter.unique_fragments_from_ring_openings), 0)
Expand Down

0 comments on commit ee27a3b

Please sign in to comment.