Molecule |
molecule<CalculationT, MoleculeProperties, AtomProperties, BondProperties, AtomTag, BondTag, MoleculeTag>
The molecule class is the basic structure of the Chemical Descriptors Library. It is parameterized on 3 fundamental components: Atom, Bond and Molecule properties. The first two properties represent the atom and bond objects. The CalculationT template parameter, is the c++ fundamental type on which calculations will be performed. AtomTag, BondTag and MoleculeTag are tags that identifies their respective properties. If you do not use the default parametes for the molecule class, you have to provide the property accessors of each property in the respective props_selectors header files that are defined in cdl/detail/
namespace morpho { namespace cdl { template < class CalculationT , class MoleculeProperties , class AtomProperties , class BondProperties , class MoleculeTag , class AtomTag , class BondTag > class molecule { public : //------------- the property types typedef AtomProperties atom_properties; typedef BondProperties bond_properties; typedef MoleculeProperties prop_type; //------------- typedef CalculationT calc_type; //------------- types contained in the property maps typedef entity<atom_properties> atom_type; typedef entity<bond_properties> bond_type; typedef entity<prop_type> molecular_props_type; //------------- tags types for the property maps typedef AtomTag atom_tag_type; typedef BondTag bond_tag_type; typedef MoleculeTag molecule_tag_type; //------------- properties of the molecular graph typedef boost::property<bond_tag_type,bond_type> EdgeProperty; // bond typedef boost::property<atom_tag_type,atom_type> VertexProperty; // atom typedef boost::property<molecule_tag_type,molecular_props_type> GraphProperty; // molecule typedef boost::adjacency_list<boost::vecS, boost::vecS, boost::undirectedS, VertexProperty, EdgeProperty, GraphProperty> graph_type; //------------- self type typedef molecule self; //------------- graph descriptors typedef graph_traits::vertex_descriptor vertex_descriptor; typedef graph_traits::edge_descriptor edge_descriptor; //------------- graph representation iterators typedef graph_traits::vertex_iterator vertex_iterator; typedef graph_traits::edge_iterator edge_iterator; typedef graph_traits::adjacency_iterator adjacency_iterator; //traverses vertices typedef graph_traits::out_edge_iterator out_edge_iterator; //traverses edges //------------- properties iterators (uses iterator adaptors from the graph repres. iterators) typedef boost::property_map_iterator_generator<atom_property_map_type, vertex_iterator>::type atom_iterator; typedef boost::property_map_iterator_generator<bond_property_map_type, edge_iterator>::type bond_iterator; typedef boost::property_map_iterator_generator<const_atom_property_map_type, vertex_iterator>::type const_atom_iterator; typedef boost::property_map_iterator_generator<const_bond_property_map_type, edge_iterator>::type const_bond_iterator; typedef boost::property_map_iterator_generator<atom_property_map_type, adjacency_iterator>::type adjacency_atom_iterator; typedef boost::property_map_iterator_generator<bond_property_map_type, out_edge_iterator>::type incident_bond_iterator; typedef boost::property_map_iterator_generator<const_atom_property_map_type, adjacency_iterator>::type adjacency_const_atom_iterator; typedef boost::property_map_iterator_generator<const_bond_property_map_type, out_edge_iterator>::type incident_const_bond_iterator; //------------- copy constructor molecule(const self& other_); //------------- constructors molecule(const graph_type& g); template <class EdgeIterator, class AtomPropIterator, class BondPropIterator> molecule(EdgeIterator ei, EdgeIterator ei_end, size_t n, AtomPropIterator api, BondPropIterator bpi, MoleculeProperties mp, bool add_implicit_h = false, bool initialize = true); molecule(const nail_juice<self>& nj, bool add_implicit_h = false, bool initialize = true); //------------- operators self& operator=(const self& rhs); //------------------------------ member functions -------------------------- // //---------- molecule // void initialize_properties(bool); const graph_type& get_graph() const; graph_type& get_graph(); // O(V+E) template <class EdgeFilterFactory, class VertexFilterFactory> graph_type get_filtered_graph(EdgeFilterFactory, VertexFilterFactory); // O(V) template <class VertexFilterFactory> graph_type get_filtered_graph(VertexFilterFactory); // O(V+E) graph_type get_filtered_graph(); // O(V+E) template <class EdgePredicate, class VertexPredicate> boost::filtered_graph<graph_type,EdgePredicate,VertexPredicate> get_boost_filtered_graph(); edges_size_type num_connections() const; // //---------- vertices // std::pair<vertex_iterator,vertex_iterator> vertices() const; std::pair<adjacency_iterator,adjacency_iterator> adjacent_vertices(vertex_descriptor) const; vertex_iterator get_vertex_iterator(const atom_iterator&) { return ai.base(); } vertex_iterator get_vertex_iterator(const const_atom_iterator&) const { return ai.base(); } // //---------- atoms // std::pair<adjacency_atom_iterator,adjacency_atom_iterator> adjacent_atoms(vertex_descriptor); std::pair<adjacency_const_atom_iterator,adjacency_const_atom_iterator> adjacent_atoms(vertex_descriptor) const; std::pair<adjacency_atom_iterator,adjacency_atom_iterator> adjacent_atoms(const atom_iterator&); std::pair<adjacency_const_atom_iterator,adjacency_const_atom_iterator> adjacent_atoms(const const_atom_iterator&) const; typename boost::call_traits<atom_type>::const_reference get_atom(vertex_descriptor) const; typename boost::call_traits<atom_type>::reference get_atom(vertex_descriptor); typename boost::call_traits<atom_type>::const_reference atom(vertices_size_type) const; typename boost::call_traits<atom_type>::reference atom(vertices_size_type); atom_iterator get_atom_iterator(vertex_iterator); const_atom_iterator get_atom_iterator(vertex_iterator) const; atom_iterator atoms_begin(); const_atom_iterator atoms_begin() const; atom_iterator atoms_end(); const_atom_iterator atoms_end() const; std::pair<atom_iterator,atom_iterator> atoms(); std::pair<const_atom_iterator,const_atom_iterator> atoms() const; vertices_size_type num_atoms() const; degree_size_type connected_atoms(vertex_descriptor) const; degree_size_type connected_atoms(const atom_iterator&); degree_size_type connected_atoms(const const_atom_iterator&) const; // //---------- edges // std::pair<edge_iterator,edge_iterator> edges() const; std::pair<out_edge_iterator, out_edge_iterator> incident_edges(vertex_descriptor) const; edge_iterator get_edge_iterator(const bond_iterator&); edge_iterator get_edge_iterator(const const_bond_iterator&) const; std::pair<edge_descriptor, bool> edge(vertex_descriptor, vertex_descriptor) const; // //---------- bonds // //gets an iterator range providing access to the incident edges of a vertex std::pair<incident_bond_iterator,incident_bond_iterator> incident_bonds(vertex_descriptor); std::pair<incident_const_bond_iterator,incident_const_bond_iterator> incident_bonds(vertex_descriptor) const; std::pair<incident_bond_iterator,incident_bond_iterator> incident_bonds(const atom_iterator&); std::pair<incident_const_bond_iterator,incident_const_bond_iterator> incident_bonds(const const_atom_iterator&) const; vertex_descriptor source(const edge_descriptor&) const; typename boost::call_traits<atom_type>::reference source(const bond_iterator&); typename boost::call_traits<atom_type>::const_reference source(const const_bond_iterator&) const; vertex_descriptor target(const edge_descriptor& ) const; typename boost::call_traits<atom_type>::reference target(const bond_iterator& ); typename boost::call_traits<atom_type>::const_reference target(const const_bond_iterator& ) const; typename boost::call_traits<bond_type>::const_reference get_bond(const edge_descriptor& ) const; typename boost::call_traits<bond_type>::reference get_bond(const edge_descriptor& ); bond_iterator get_bond_iterator(edge_iterator ei); const_bond_iterator get_bond_iterator(edge_iterator ei) const; bond_iterator bonds_begin(); const_bond_iterator bonds_begin() const; bond_iterator bonds_end(); const_bond_iterator bonds_end() const; std::pair<bond_iterator,bond_iterator> bonds(); std::pair<const_bond_iterator,const_bond_iterator> bonds() const; boost::compressed_pair<bond_type&, bool> bond(vertex_descriptor , vertex_descriptor ); boost::compressed_pair<const bond_type&, bool> bond(vertex_descriptor , vertex_descriptor ) const; boost::compressed_pair<bond_type&, bool> bond(const atom_iterator& , const atom_iterator& ); boost::compressed_pair<const bond_type&, bool> bond(const const_atom_iterator& , const const_atom_iterator& ) const; unsigned num_bonds() const; unsigned num_bonds(vertex_descriptor ) const; unsigned num_bonds(const atom_iterator& ); unsigned num_bonds(const const_atom_iterator& ) const; // //---------- structure modification // edge_descriptor add_atom(vertex_descriptor, const char* symbol); edge_descriptor add_atom(vertex_descriptor, std::string symbol); edge_descriptor add_atom(atom_iterator, const char*); edge_descriptor add_atom(atom_iterator, std::string); vertex_descriptor add_atom(const char*); void clear_atom(vertex_descriptor); void remove_atom(vertex_descriptor); template <class Predicate> void remove_atoms_if(Predicate, bool); void replace_atom(vertex_descriptor, atom_iterator); void replace_atom(atom_iterator, atom_iterator); void replace_atom(vertex_descriptor, const char* ); void replace_atom(atom_iterator, const char* ); std::pair<edge_descriptor, bool> add_bond(vertex_descriptor, vertex_descriptor); std::pair<edge_descriptor, bool> add_bond(atom_iterator, atom_iterator); void remove_bond(vertex_descriptor, vertex_descriptor); void remove_bond(atom_iterator, atom_iterator); void remove_bond(edge_descriptor); void remove_bond(bond_iterator); templatevoid remove_incident_bonds_if(vertex_descriptor, UnaryPredicate); template void remove_bonds_if(UnaryPredicate); //------------------------------ member variables -------------------------- private : graph_type m_molecular_graph; }; // class molecule } // namespace cdl } // namespace morpho
Parameter | Description | Default | Model of / Requirements |
---|---|---|---|
CalculationT | Type to perform the calculations | double | fundamental type |
AtomProperties | Object that contains the atomic properties | boost::property< atom_static_propsS, ::morpho::atomic_props::atomic_properties<>const*,boost::property< atom_rw_propsS, rw_atomic_properties<> > > | boost::property<> |
BondProperties | Object that contains the bond properties | boost::property< bond_propsS, bond_properties<> > | boost::property<> |
MoleculeProperties | Objects that contains the molecul props | boost::property< mol_propsS, molecular_properties<> > | boost::property<> |
AtomTag | Tag for the atomic props | atom_t | constructible |
BondTag | Tag for the bonds props | bond_t | constructible |
MoleculeTag | Tag for the molecule props | mol_t | constructible |
All three properties objects (atom, bond and molecule) must model boost::property<> (defined in boost/pending/property.hpp). Boost property class is just an interface to access properties in a generic way. It is defined as follow:
template <class Tag, class T, class Base = no_property>
struct property : public Base
It has 3 template parameters: Tag: tag of the property, T: the property by itself,
and Base: can be other boost::property<>. With the Base parameter, it is possible to concatenate an
arbitrary number of properties, and have the same accessing interface.
Bond and molecule properties can be of arbitrary type, and can have arbitrary tags. However,
the AtomicProperties template parameter has the following requierement:
Must have at least 1 property, one of these must have the tag atom_static_propsS and must be a pointer to a constant morpho::atomic_props::atomic_properties<> class. Here CDL retrieves the atomic properties that are constant. |
atom_properties
The type of your Vertex property: in CDL context, properties for an atom
bond_properties
The type of your Edge property: in CDL context, properties for a bond
prop_type
The type of the property for the molecule itself
calc_type
Fundamental type used for calculations on the molecule
atom_type
The type of entity<AtomProperties>. This is the type atom_iterator points to.
bond_type
The type of entity<BondProperties>. This is the type bond_iterator points to.
molecular_props_type
The type of entity<MoleculeProperties>. The properties of the molecule.
atom_tag_type
A tag type for the atom properties
bond_tag_type
A tag type for the bond properties
molecule_tag_type
A tag type for the molecule properties
graph_type
The type of the graph member variable
vertex_descriptor
The type for the vertex descriptors associated with the adjacency_list. A vertex descriptor is the type of the element that a vertex iterator points to
edge_descriptor
The type for the edge descriptors associated with the adjacency_list. An edge descriptor is the type of the element that an edge iterator points to
vertex_iterator
The type for the iterators returned by vertices()
edge_iterator
The type for the iterators returned by edges()
out_edge_iterator
The type for the iterators returned by out_edges()
adjacency_iterator
The type for the iterator returned by adjacent_vertices(). The adjacency_iterator models the same iterator concept as out_edge_iterator
atom_iterator
The type for the iterator returned by atoms_begin(), atoms_end() and atoms()
bond_iterator
The type for the iterator returned by bonds_begin(), bonds_end() and bonds()
const_atom_iterator
Same as the atom_iterator, but points to a constant atom
const_bond_iterator
Same as the bond_iterator, but points to a constant bond
adjacency_atom_iterator
The type for the iterator returned by adjacent_atoms()
incident_bond_iterator
The type for the iterator returned by incident_bonds()
adjacency_const_atom_iterator
Same as the adjacency_atom_iterator, but points to a constant atom
incident_const_bond_iterator
Same as the incident_bond_iterator, but points to a constant bond
molecule(const molecule&)
Copy constructor. Creates copies of all internal properties, including atoms and bonds. Also copies the the molecule property.
molecule(const graph_type&)
Creates a molecule out of a graph_type. This copies all the internal properties of the graph, and initializes the molecule properties.
template <class EdgeIterator, class AtomPropIterator, class BondPropIterator> molecule(EdgeIterator ei, EdgeIterator ei_end, size_t n, AtomPropIterator api, BondPropIterator bpi, MoleculeProperties mp, bool add_implicit_h = false, bool initialize = true)
Creates a molecule object with n atoms, and bonds indicated in the bond's list given by the iterator range [ei, ei_end). EdgeIterator value is of type std::pair<size_t, size_t>, where size_t corresponds to atom indices, which must fall in the range [0, n): is the connectivity table. AtomPropIterator and BondPropIterator model InputIterator, and point to constructed atom and bond property objects respectively. MoleculeProperties is a constructed molecule properties object. Flag initialize tells wheter to initialize all internal properties or not.
molecule(const nail_juice<self>& nj, bool add_implicit_h = false, bool initialize = true)
Creates a molecule object and attaches the information that comes with the nail_juice to the properties if the newly constructed molecule. Please refer to nail_juice documentation when using this constructor. Flag initialize tells wheter to initialize all internal properties or not.
self& operator=(const self& rhs)
Copies rhs molecular graph into the this molecular graph member. Returns *this.
void initialize_properties(bool add_implicit_h)
Initializes all internal properties to the calculated values. It does, between other calculations, the assignment of the sssr verties and edges, assignment of flags, calculation of distance matrix, etc.
const graph_type& get_graph() constand
graph_type& get_graph()
Returns the graph object
template <class EdgeFilterFactory, class VertexFilterFactory> graph_type get_filtered_graph(EdgeFilterFactory, VertexFilterFactory)
Returns a filtered graph.
Complexity: exactly V applications
of the unary predicate VertexFilter plus E applications of
unary predicate EdgeFilter
Requirenment for EdgeFilterFactory and VertexFilterFactory:
have a member function called make_filter() which accepts a molecule, and returns a filter. A filter is a unary predicate that accepts an atom or a bond (depending for which objects the filter is) that discriminates which atom-bond should be filtered. |
template <class VertexFilterFactory> graph_type get_filtered_graph(VertexFilterFactory)Returns an atom-filtered graph. It preserves all bonds
graph_type get_filtered_graph()Returns a graph without H-bonds and without H-atoms
graph_type get_filtered_graph()Returns a pure boost::filtered_graph.
size_t num_connections() const
Returns the number of connections the molecule has. Note that double and triple bonds are counted as just one connection. If you want the total number of bonds, use num_bonds() instead
std::pair<vertex_iterator,vertex_iterator> vertices() constReturns an iterator-range providing access to the vertex set of graph the molecular graph
std::pair<adjacency_iterator,adjacency_iterator> adjacent_vertices(vertex_descriptor v) constReturns an iterator-range providing access to the adjacent vertices of vertex v
vertex_iterator get_vertex_iterator(const atom_iterator& ai) vertex_iterator get_vertex_iterator(const const_atom_iterator& ai) constReturns the vertex_iterator for the vertex from which ai is mapped
std::pair<adjacency_atom_iterator,adjacency_atom_iterator> adjacent_atoms(vertex_descriptor v) std::pair<adjacency_const_atom_iterator,adjacency_const_atom_iterator> adjacent_atoms(vertex_descriptor v) const
Returns an iterator range providing access to the adjacent atoms of the atom mapped from v
std::pair<adjacency_atom_iterator,adjacency_atom_iterator> adjacent_atoms(const atom_iterator& ai) std::pair<adjacency_const_atom_iterator,adjacency_const_atom_iterator> adjacent_atoms(const const_atom_iterator& ai) const
Returns an iterator range providing access to the adjacent atoms of the atom pointed to by ai
boost::call_traits<atom_type>::const_reference get_atom(vertex_descriptor v) const boost::call_traits<atom_type>::reference get_atom(vertex_descriptor v)Returns the atom that is mapped from vertex v
boost::call_traits<atom_type>::const_reference atom(vertices_size_type n) const boost::call_traits<atom_type>::reference atom(vertices_size_type n)Returns the 'n'th atom in the graph's vertex list
atom_iterator get_atom_iterator(vertex_iterator vi) const_atom_iterator get_atom_iterator(vertex_iterator vi) constConverts a vertex_iterator to an atom_iterator
atom_iterator atoms_begin() const_atom_iterator atoms_begin() constReturns an atom_iterator that points to the begin of the atom sequence
atom_iterator atoms_end() const_atom_iterator atoms_end() constReturns an atom_iterator that points to the end of the atom sequence
std::pair<atom_iterator,atom_iterator> atoms() { std::pair<const_atom_iterator,const_atom_iterator> atoms() const {Returns an iterator range providing access to the atoms of the molecule
vertices_size_type num_atoms() constReturns the number of atoms in the molecule
size_t connected_atoms(vertex_descriptor u) constReturns the number of vertices connected to vertex u
degree_size_type connected_atoms(const atom_iterator& ai) degree_size_type connected_atoms(const const_atom_iterator& ai) const {Returns the number of atoms connected to the atom pointed to by ai
std::pair<edge_iterator,edge_iterator> edges() constReturns an iterator range providing access to the edges set of the molecular graph
std::pair<out_edge_iterator, out_edge_iterator> incident_edges(vertex_descriptor v) constReturns an iterator range providing access to the incident edges of vertex v
edge_iterator get_edge_iterator(const bond_iterator& bi) edge_iterator get_edge_iterator(const const_bond_iterator& bi) constConverts a bond iterator to an edge_iterator
std::pair<edge_descriptor, bool> edge(vertex_descriptor u, vertex_descriptor v) constReturns an edge connecting vertex u to vertex v
std::pair<incident_bond_iterator,incident_bond_iterator> incident_bonds(vertex_descriptor v) std::pair<incident_const_bond_iterator,incident_const_bond_iterator> incident_bonds(vertex_descriptor v) constReturns an iterator range providing access to the incident bonds of the atom mapped from vertex v
std::pair<incident_bond_iterator,incident_bond_iterator> incident_bonds(const atom_iterator& ai) std::pair<incident_const_bond_iterator,incident_const_bond_iterator> incident_bonds(const const_atom_iterator& ai) constReturns an iterator range providing access to the incident bonds of the atom pointed to by ai
vertex_descriptor source(const edge_descriptor& e) constReturns the vertex_descriptor for the source of the edge_descriptor e
boost::call_traits<atom_type>::reference source(const bond_iterator& bi) {Returns a reference to the atom that is the source of the bond pointed to by bi
boost::call_traits<atom_type>::const_reference source(const const_bond_iterator& bi) const {
Returns a const reference to the atom that is the source of the bond pointed to by bi
vertex_descriptor target(const edge_descriptor& e) const {Returns the vertex_descriptor for the target of the edge_descriptor e
boost::call_traits<atom_type>::reference target(const bond_iterator& bi) {Returns a reference to the atom that is the target of the bond pointed to by bi
boost::call_traits<atom_type>::const_reference target(const const_bond_iterator& bi) const {Returns a const reference to the atom that is the target of the bond pointed to by bi
boost::call_traits<bond_type>::reference get_bond(const edge_descriptor& e) {Returns a reference to the bond object that is mapped from the e edge
boost::call_traits<bond_type>::const_reference get_bond(const edge_descriptor& e) constReturns a const reference to the bond object that is mapped from the e edge
bond_iterator get_bond_iterator(edge_iterator ei) const_bond_iterator get_bond_iterator(edge_iterator ei) constReturns a bond_iterator which points to a bond (and a constant bond respectively) that is mapped from the edge pointed to by ei
bond_iterator bonds_begin() const_bond_iterator bonds_begin() constReturns an bond_iterator that points to the begin of the atom sequence
bond_iterator bonds_end() const_bond_iterator bonds_end() constReturns an bond_iterator that points to the end of the atom sequence
std::pair<bond_iterator,bond_iterator> bonds() std::pair<const_bond_iterator,const_bond_iterator> bonds() constReturns an iterator range providing access to the bonds of the molecule
boost::compressed_pair<bond_type&, bool> bond(vertex_descriptor v1, vertex_descriptor v2) boost::compressed_pair<const bond_type&, bool> bond(vertex_descriptor v1, vertex_descriptor v2) constReturns a reference to a bond that is mapped from the edge that creates the connection between v1 and v2. If the bond doesn't exist, it returns the first bond, but with the second parameter of the pair set to false
boost::compressed_pair<bond_type&, bool> bond(const atom_iterator& ai1, const atom_iterator& ai2) boost::compressed_pair<const bond_type&, bool> bond(const const_atom_iterator& ai1, const const_atom_iterator& ai2) constReturns a reference to a bond that is mapped from the edge that creates the connection between the atoms pointed to by ai1 and ai2. If the bond doesn't exist, it returns the first bond, but with the second parameter of the pair set to false
unsigned num_bonds() constReturns the number of bonds in the molecule. This function counts double and triple bonds
unsigned num_bonds(vertex_descriptor v) constReturns the number of bonds that connects to the atom that is mapped from vertex v
unsigned num_bonds(const atom_iterator& ai) unsigned num_bonds(const const_atom_iterator& ai) constReturns the number of bonds that connects to the atom that is pointed by ai
edge_descriptor add_atom(vertex_descriptor v, const char* symbol) edge_descriptor add_atom(vertex_descriptor v, std::string symbol)Adds an atom of symbol symbol to the molecule, and connects it to vertex v. The function returns the edge descriptor for that connection. This function invalidates out_edge_iterator and incident_bond_iterator for vertex v
edge_descriptor add_atom(atom_iterator ai, const char* symbol) edge_descriptor add_atom(atom_iterator ai, std::string symbol)Adds an atom of symbol symbol to the molecule, and connects it to the atom pointed to by ai. The function returns the edge descriptor for that connection. The function invalidates out_edge_iterator and incident_bond_iterator for vertex v
vertex_descriptor add_atom(const char* symbol)Adds a single atom to the molecule. The function doesn't connect the atom to another atom. The function returns the index for the new vertex.
void clear_atom(vertex_descriptor v);Removes all incident edges of vertex v. The atom still appears in the atom's set. Causes the same on iterators and descriptors as the remove_edge() function.
void remove_atom(vertex_descriptor v);Removes all incident edges of vertex v, then removes v itself. This function invalidates all descriptors and iterators for the graph.
template <class Predicate> void remove_atoms_if(Predicate p, bool reinitialize_props = true);Removes all incident edges and then removes the atoms that fulfill the predicate. This function invalidates all descriptors and iterators for the graph. The bool flag tells the algorithm if the properties should be reinitialized. class Predicate accepts an atom.
void replace_atom(vertex_descriptor v, atom_iterator other_atom)Replaces atom mapped from vertex v with an atom pointed to by other_atom.
void replace_atom(atom_iterator this_atom, atom_iterator other_atom)Replaces atom pointed to by this_atom, with the atom pointed to by other_atom.
void replace_atom(vertex_descriptor v, const char* new_atom_symbol)Replaces atom mapped from vertex v with a new atom of symbol new_atom_symbol.
void replace_atom(atom_iterator ai, const char* new_atom_symbol)Replaces atom pointed to by ai with a new atom of symbol new_atom_symbol.
std::pair<edge_descriptor, bool> add_bond(vertex_descriptor u, vertex_descriptor v)Adds edge (u,v) to the molecular graph. This operation will invalidate any out_edge_iterator for vertex u and v, and of course, the incident_bond_iterator for the atoms mapped from vertex u and v. bool is a flag that indicates if the the bond could be added. The functions returns the edge descriptor for the added edge. If the flag is false, the returned edge descriptor points to the already existing edge.
std::pair<edge_descriptor, bool> add_bond(atom_iterator ai1, atom_iterator ai2)Adds a bond that connects atoms pointed to by ai1 and ai2. This operation will invalidate out_edge_iterators and incident_bond_iterators for the atom.
void remove_bond(vertex_descriptor u, vertex_descriptor v)Removes edge (u,v) from the molecular graph. The operation causes edge descriptors or iterators pointing to the incident edges or bonds of vertex v and u to become invalid. Also invalidates iterators that point to the edge list of vertices u and v.
void remove_bond(atom_iterator ai1, atom_iterator ai2)Removes bond that connects atoms pointed to by ai1 and ai2. The operation invalidates all iterators that points to the bond or edge.
void remove_bond(edge_descriptor e)Removes edge e from the molecular graph. Invalidates all iterators that points to that edge.
void remove_bond(bond_iterator bi)Removes bond that is pointed to by bi. Invalidates all iterators that points to the bond or edge.
templateRemoves the incident bonds of vertex v if it satisfies predicate pred. Predicate accepts an edge descriptor.void remove_incident_bonds_if(vertex_descriptor v, UnaryPredicate pred)
templateRemoves all bonds from the molecular graph if satisfies the unary predicate pred. Unary predicate accepts an edge descriptor.void remove_bonds_if(UnaryPredicate pred)