Internal types and functions
See also the documentations of PeriodicGraphs.jl
and of PeriodicGraphEmbeddings.jl
Types
CrystalNets.Crystal
— TypeCrystal
Intermediate representation of a crystal, retaining information on the cell, and the fractional placement of the atoms and their type, as well as the residues which will be used as vertices for the computation of the underlying topology.
CrystalNets.Clusters
— TypeClusters
Classification of the atoms of a crystalline framework in different clusters. For simple crystals, every atom is its own cluster. For a MOF, a cluster is a SBU, which can be either organic or inorganic.
CrystalNets.CollisionNode
— TypeCollisionNode
Store the structure of a collision node through the subgraph g
extracted with only the edges bond to the vertices in the node.
The num
field corresponds to the number of vertices in g
that collide in the initial graph, while the neighs
field stores the indices of their neighbors out of the collision site.
The colliding vertices are the num
first vertices of g
, the others are the neighbors. In particular, nv(g) == num + length(neighs)
and g[(num+1):(nv(g))]
has no edge.
CrystalNets.CIF
— TypeCIF
Representation of a .cif file.
Core topology functions
CrystalNets.topological_key
— Functiontopological_key(net::CrystalNet)
Return a unique topological key for the net, which is a topological invariant of the net (i.e. it does not depend on its initial representation).
CrystalNets.CRYSTALNETS_ARCHIVE
— Constantconst CRYSTALNETS_ARCHIVE::Dict{String,String}
The archive used to recognize known topologies.
You probably don't need to access it directly: rely on recognize_topology
to read and the various archive functions like add_to_current_archive!
to write.
CrystalNets.minimize
— Functionminimize(net::CrystalNet, [collisions::Vector{CollisionNode}])
Return a CrystalNet representing the same net as the input, but in a unit cell. If collisions
is given, also return the corresponding collisions after minimization.
The computed unit cell may depend on the representation of the input, i.e. it is not topologicallly invariant.
minimize(net::Crystal{Nothing})
Return a Crystal representing the same crystal as the input, but in a unit cell.
CrystalNets.candidate_key
— Functioncandidate_key(net::CrystalNet, u, basis, minimal_edgs)
Given the net, a candidate u => basis
where u
is the origin and basis
the triplet of axes, and minimal_edgs
the last minimal key (for the pseudo-lexicographical order used), extract the key corresponding to the current candidate.
The key is the lexicographically ordered list of edges of the graph when its vertices are numbered according to the candidate. The ordering of keys first compares the list of edges disregarding the offsets, and then only compares the offsets if the rest is identical.
If the key is larger or equal to minimal_edgs
, early stop and return two empty lists. Otherwise, the extracted key is the current best: return the vmap between the initial vertices and their ordered image in the candidate, as well as the key.
See also: find_candidates
CrystalNets.possible_translations
— Functionpossible_translations(c::Union{CrystalNet,CrystalNet})
Return a list of tuples (nz, i_max_den, max_den, t)
where
t
is a translation mapping at the origin vertex to another one in the unit cell.max_den
is the maximum denominator in theD
coefficients oft
.i_max_den
is the index.nz
is the number of zeros int
.
The list is guaranteed to contain all the possible valid translations but may contain some invalid translations.
See also: find_all_valid_translations
, PeriodicGraphEmbeddings.check_valid_symmetry
CrystalNets.find_all_valid_translations
— Functionfind_all_valid_translations(c::Union{Crystal,CrystalNet{D}}, collisions) where D
Return a D
-tuple of list of tuples (i_max_den, max_den, t)
(see possible_translations
for interpretation) where the n
-th list contains all valid translations of the net having exactly n-1
zeros.
A translation is valid if it maps exactly each vertex to a vertex and each edge to an edge.
See also: possible_translations
, PeriodicGraphEmbeddings.check_valid_symmetry
CrystalNets.minimal_volume_matrix
— Functionminimal_volume_matrix(translations::NTuple{D}) where {D}
Given the output of find_all_valid_translations
, compute the transformation that allows reducing the net to its minimal cell.
CrystalNets.reduce_with_matrix
— Functionreduce_with_matrix(c::CrystalNet, mat, collisions)
Given the net and the output of minimal_volume_matrix
computed on the valid translations of the net, return the new net representing the initial net in the computed unit cell.
CrystalNets.partition_by_coordination_sequence
— Functionpartition_by_coordination_sequence(graph, symmetries::AbstractSymmetryGroup=NoSymmetryGroup(graph))
Partition the vertices of the graph into disjoint categories, one for each coordination sequence. The partition is then sorted by order of coordination sequence. This partition does not depend on the representation of the graph. The optional argument vmaps
is a set of permutations of the vertices that leave the graph unchanged. In other words, vmaps
is a set of symmetry operations of the graph.
Return the categories and a list of unique representative for each symmetry class.
CrystalNets.find_candidates
— Functionfind_candidates(net::CrystalNet{D}, collisions::Vector{CollisionNode}) where D
Return a non-empty set of candidates u => basis
where u
is a vertex and basis
is matrix whose columns are D
linearly independent euclidean embeddings of edges. The returned set is independent of the representation of the graph used in net
.
Also return a category_map
linking each vertex to its category number, as defined by partition_by_coordination_sequence
See also: candidate_key
CrystalNets.extract_through_symmetry
— Functionextract_through_symmetry(candidates::Dict{Int,Vector{SMatrix{3,3,T,9}}}, symmetries::AbstractSymmetryGroup) where T
Given the candidates and the list of symmetries of the net, return the flattened list of candidates after removing candidates that are symmetric images of the kept ones.
CrystalNets.find_initial_candidates
— Functionfind_initial_candidates(net::CrystalNet{D}, candidates_v, category_map) where D
Given the net, a list of vertices in a given category and the category_map
, return a list of pairs u => (basis, cats)
where u ∈ candidates_v
, basis
is a D
-rank matrix made by juxtaposing the euclidean embeddings of outgoing edges from u
, and cats
are the categories of the respective neighbors of u
.
If the basis
corresponding to vertex u
is not of rank D
, it is not included in the returned list (for instance, if all outgoing edges of a vertex are coplanar with D == 3
).
CrystalNets.find_candidates_onlyneighbors
— Functionfind_candidates_onlyneighbors(net::CrystalNet{D}, candidates_v, category_map) where D
Given the net, a list of vertices in a given category and the category_map
, return a Dict whose pairs u => matv
are such that u ∈ candidates_v
and matv
is a list of unique invertible matrices of size D
whose columns are euclidean embeddings of outgoing edges from u
. Each such matrix has a category, defined by the D
-uplet of categories of each corresponding outneighbor of u
: the returned Dict is such that all the matrices belonging to all matv
share the same category.
The returned Dict is empty iff find_initial_candidates(net, candidates_v, category_map)
is empty.
CrystalNets.find_candidates_fallback
— Functionfind_candidates_fallback(net::CrystalNet3D, reprs, othercats, category_map)
Return candidates in the same form as find_candidates_onlyneighbors
except that only two edges start from u
and one does not.
Input
CrystalNets.parse_cif
— Functionparse_cif(file_path)
Parse a CIF file and return a dictionary where each identifier (without the starting '_' character) is linked to its value. Values are either a string or a vector of string (if defined in a loop).
CrystalNets.CIF
— MethodCIF(file_path::AbstractString)
Make a CIF object out of the parsed file.
CrystalNets.check_collision
— Functioncheck_collision(pos, mat)
Given a list of fractional coordinates pos
and the matrix of the unit cell mat
, return a list of atoms that are suspiciously close to another atom of the list. For each collision site, only one atom is not present in the returned list.
CrystalNets.fix_atom_on_a_bond!
— Functionfix_atom_on_a_bond!(graph, pos, mat)
Remove bonds that are intercepted by an atom.
CrystalNets.least_plausible_neighbours
— Functionleast_plausible_neighbours(Δs, n)
Find the positions of the n least probable neighbours of an atom, given the list Δs of the distance between their position and that of the atom.
This function is highly empirical and should not be considered utterly reliable.
CrystalNets.fix_valence!
— Functionfix_valence!(graph::PeriodicGraph3D, pos, types, passH, passO, passCN, mat, ::Val{dofix}, options) where {dofix}
Attempt to ensure that the coordinence of certain atoms are at least plausible by removing some edges from the graph. These atoms are H, halogens, O, N and C. if dofix
is set, actually modify the graph; otherwise, only emit a warning. In both cases, return a list of atoms with invalid coordinence.
CrystalNets.sanitize_removeatoms!
— Functionsanitize_removeatoms!(graph::PeriodicGraph3D, pos, types, mat, options)
Special heuristics to remove atoms that seem to arise from an improper cleaning of the file. Currently implemented:
- C atoms suspiciously close to metallic atoms.
- One of two identical metallic atoms suspiciously close to one another
TODO:
- O atoms with 4 coplanar bonds (warning only).
CrystalNets.remove_triangles!
— Functionremove_triangles!(graph::PeriodicGraph3D, pos, types, mat, toinvestigate=collect(edges(graph)))
In a configuration where atoms A, B and C are pairwise bonded, remove the longest of the three bonds if it is suspicious (too large and too close to the third atom).
CrystalNets.remove_homoatomic_bonds!
— Functionremove_homoatomic_bonds!(graph::PeriodicGraph, types, targets, reduce_homometallic)
Remove from the graph all bonds of the form X-X where X is an atom in targets
.
Also remove all such bonds where X is a metal if the two bonded atoms up to third neighbours otherwise, and if reduce_homometallic
is true
.
CrystalNets.sanity_checks!
— Functionsanity_checks!(graph, pos, types, mat, options)
Perform some sanity checks to ensure that the detected bonds are not obviously wrong because they are either too short or too long.
Crystal and CIF handling
CrystalNets.remove_partial_occupancy
— Functionremove_partial_occupancy(::CIF)
Only keep one atom per atom site.
CrystalNets.prune_collisions
— Functionprune_collisions(::CIF)
For each site where there are atoms suspiciously close to one another, remove all but one of them. This arises for example when all the possible positions of at atom are superposed in the CIF file, typically for a solvent which should be disregarded anyway.
CrystalNets.expand_symmetry
— Functionexpand_symmetry(::CIF)
Applies all the symmetry operations listed in the CIF file to the atoms and the bonds.
CrystalNets.trim_monovalent
— Functiontrim_monovalent(crystal)
Repeatedly remove monovalent atoms from the crystal until none is left.
CrystalNets.trimmed_crystal
— Functiontrimmed_crystal(c::Crystal{Nothing})
Rebuild the crystal after trimming its graph according to trim_topology
.
CrystalNets.trim_topology
— Functiontrim_topology(graph::PeriodicGraph)
Return a pair (vmap, newgraph)
extracted from the input by removing vertices of valence lower or equal to 1, and by replacing vertices of valence 2 by edges, until convergence. The only exceptions are vertices only bonded to their representatives of another cell: those will not be replaced by edges even if their valence is 2, since this latter case indicates an irreducible trivial 1-dimensional topology.
vmap
maps the vertices of newgraph
to their counterpart in graph
.
Bond guessing
CrystalNets.guess_bonds
— Functionguess_bonds(pos, types, mat, options)
Return the bonds guessed from the positions, types and cell matrix, given as a Vector{Vector{Tuple{Int,Float32}}}
.
The i
-th entry of the list is a list, whose entries are of the form (j, dist)
which indicates that the representatives of vertices i
and j
distant of at most dist
are bonded together.
CrystalNets.edges_from_bonds
— Functionedges_from_bonds(bonds::Vector{Vector{Tuple{Int,Float32}}}, mat, pos)
Given a bond list bonds
containing triplets (a, b, dist)
where atoms a
and b
are bonded if their distance is lower than dist
, the 3×3 matrix of the cell mat
and the Vector{SVector{3,Float64}} pos
whose elements are the fractional positions of the atoms, extract the list of PeriodicEdge3D corresponding to the bonds. Since the adjacency matrix wraps bonds across the boundaries of the cell, the edges are extracted so that the closest representatives are chosen to form bonds.
Clustering algorithm
CrystalNets.find_sbus!
— Functionfind_sbus!(crystal::Crystal, kinds::ClusterKinds=default_sbus)
Recognize SBUs using heuristics based on the atom types corresponding to the AllNodes
clustering algorithm.
CrystalNets.regroup_sbus
— Functionregroup_sbus(graph::PeriodicGraph3D, classes::AbstractVector{<:Integer},
isolate=Int[])
Given a classification of vertices into classes, separate the vertices into clusters of contiguous vertices belonging to the same class.
isolate
is a list where each atom is separated from the rest of its class. Once all such atoms of its class are isolated, we look for the connected components of non-isolated atoms in that class. If such a component has only one neighbours which is an isolated atom, it is added to the vertex of the isolated atom.
CrystalNets.regroup_paddlewheel!
— Functionregroup_paddlewheel!(graph, clusters::Clusters, types, periodicsbus)
Identify paddle-wheel patterns made of two opposite SBUs and regroup them into one.
CrystalNets.split_sbu!
— Functionsplit_sbu!(sbus, graph, i_sbu, classes)
Split SBU number i_sbu
into new SBUs according to the updated classes
. The first argument sbus
is modified in-place. Return the list of newly-created periodic SBUs, if any.
CrystalNets.reclassify!
— Functionreclassify!(sbus, newperiodicsbus, newclass, graph, types, classof, i_sbu)
Reclassify the atoms of sbus.sbus[i_sbu])
according to the following algorithm:
- Let's call "target atom" any atom of type
typ
whereclassof[typ] == deg
and eitherdeg == 0
ordeg > 0
and the degree of the atom isdeg
. - Assign a new SBU for each target atom (one new per atom).
- Look at the connected components of atoms in the SBU which are not target atoms. For each connected component that is finite (0-dimensional) and has only one neighbor which is a target atom, put that component in the same SBU as the neighbor.
CrystalNets.add_to_newclass!
— Functionadd_to_newclass!(classes, graph, sbus, new_class, v, types, noneighborof)
Set the class of v
to new_class
. Then, grow the newly created class by adding connected components of the SBU of v
such that the new class does not become periodic and does not contain any vertex that is a neighbor of a vertex whose type is in noneighborof
.
If types === nothing
, disregard the condition on noneighborof
.
CrystalNets.group_cycle
— Functiongroup_cycle(organiccycle, types, graph)
Return a list of Vector{PeriodicVertex3D} where each sublist consists in atoms belonging to the same cycle, and which should thus belong to the same vertex eventually.
CrystalNets.collapse_clusters
— Functioncollapse_clusters(crystal::Crystal)
Return the list of crystals corresponding to the input where each cluster has been transformed into a new vertex, for each targeted clustering.
CrystalNets.pem_to_pe
— Functionpem_to_pe(cryst::Crystal{Nothing})
Convert PEM
result to PE
by removing all metallic sbus.
CrystalNets.allnodes_to_singlenodes
— Functionallnodes_to_singlenodes(cryst::Crystal)
Convert AllNodes
result to SingleNodes
by collapsing all points of extension clusters bonded together into a new organic cluster.
Unstable nets
CrystalNets.shrink_collisions
— Functionshrink_collisions(net::CrystalNet, collisions)
Remove all colliding vertices and replace them by one new vertex per CollisionNode
, whose neighbours are that of the vertices within.
CrystalNets.order_collision
— Functionorder_collision(graph::PeriodicGraph, colliding)
Given collision nodes (in the form of the corresponding list of colliding vertices), find an ordering of them which is independent of the current ordering of these vertices and of vertices which are neither in the collision node nor any of its neighbours. Return a this ordering and a priority list for each colliding vertex, or two empty lists if it fails.
This function assumes that no vertex in the node has a neighbour in another collision node and that there are no two representatives of the same vertex that are neighbour to some vertices in the range.
CrystalNets.expand_collisions
— Functionexpand_collisions(collisions::Vector{CollisionNode}, graph::PeriodicGraph, vmap)
Expand each collision node into the appropriate number of vertices so that the resulting graph is isomorphic to the initial one, in a manner that only depends on the current graph. Return the resulting graph.
vmap
is the map of vertices between initial_graph
(with collapsed collision nodes) and graph
Also return the permutations of nodes of graph
prior to expansion.
CrystalNets.collision_nodes
— Functioncollision_nodes(c::CrystalNet)
Check that the net is stable, i.e. that no two vertices have the same equilibrium placement.
A net is still considered stable if the collisions in equilibrium placement cannot lead to different topological genomes. In practice, this happens when: A) there is no edge between two collision sites and B) there is no edge between a collision site and two representatives of the same vertex and C) for each collision site, the site is made of at most 4 vertices
In this case, return the CollisionNodeList
with the corresponding CollisionNode
s, the list being empty if the net is truly stable. Otherwise, return nothing
.
Also return an updated net where the vertices in a CollisionNode
are collapsed into a new vertex, appearing after the non-colliding vertices.
Archives
CrystalNets.make_archive
— Functionmake_archive(path, destination=nothing, verbose=false)
Make an archive from the files located in the directory given by path
and export it to destination
, if specified. Each file of the directory should correspond to a unique topology: if a topology is encountered multiple times, it will be assigned the name of the latest file that bore it.
The archive can then be used with change_current_archive!(destination; validate=false)
for instance.
Utils
CrystalNets.@toggleassert
— Macro@toggleassert expression
Internal macro used to assert and expression conditionally on a build-time constant. To toggle on or off these assertions, the constant has to be modified in the source code and the module rebuilt afterwards.
CrystalNets.check_dimensionality
— Functioncheck_dimensionality(c::CrystalNet)
Check that the dimensionality of the net (i.e. the number of independent axes along which it is periodic) is equal to D
, or throw a DimensionMismatch otherwise.
Other
CrystalNets.guess_topology
— Functionguess_topology(path, options::Options)
guess_topology(path; kwargs...)
Tries to determine the topology of the file at path
by passing various options (starting from the provided options if any) until finding a known topology. If none is found, return the topological genome encountered most often through the variation of options.
CrystalNets.guess_topology_dataset
— Functionguess_topology_dataset(path, save, autoclean, showprogress, options::Options)
guess_topology_dataset(path; save=true, autoclean=true, showprogress=true, kwargs...)
Given a path to a directory containing structure input files, guess the topology of each structure within the directory using guess_topology
. Return a dictionary linking each file name to the result. The result is the corresponding topology name, if known, or the topological genome preceded by an "UNKNOWN" mention otherwise. In case of error, the result is the exception preceded by a "FAILED with" mention. Finally, if the input does not represent a periodic structure, the result is "0-dimensional".
It is strongly recommended to toggle warnings off (through toggle_warning
) and not to export any file since those actions may critically reduce performance, especially for numerous files.
The save
and autoclean
arguments work identically to their counterpart for determine_topology_dataset
.
If showprogress
is set, a progress bar will be displayed representing the number of processed files.
CrystalNets.recognize_topology
— Functionrecognize_topology(g::PeriodicGraph, arc=CRYSTALNETS_ARCHIVE)
recognize_topology(genome::AbstractString, arc=CRYSTALNETS_ARCHIVE)
Attempt to recognize a topological genome from an archive of known genomes.
This function does a simple lookup, not any kind of topology computation. To identify the topology of a PeriodicGraph
or a CrystalNet
x
, query topological_genome(x)
instead.
CrystalNets.total_interpenetration
— Functiontotal_interpenetration(itr::InterpenetratedTopologyResult, clustering::Union{Nothing,_Clustering}=nothing)
Return a Dict{TopologicalGenome,Int}
that links each topology to the number of interpenetrated nets having that topology (catenation included) for the given clustering
. If clustering
is nothing
, all possible topologies will be studied.
Example
See InterpenetratedTopologyResult
for reference on these examples.
julia> g = PeriodicGraph("2 1 1 0 2 2 2 0 1 2 2 1 0");
julia> topologies1 = topological_genome(g)
2 interpenetrated substructures:
⋅ Subnet 1 → (2-fold) UNKNOWN 1 1 1 1
⋅ Subnet 2 → sql
julia> CrystalNets.total_interpenetration(topologies1)
Dict{TopologicalGenome, Int64} with 2 entries:
UNKNOWN 1 1 1 1 => 2
sql => 1
julia> mof14 = joinpath(dirname(dirname(pathof(CrystalNets))), "test", "cif", "MOFs", "MOF-14.cif");
julia> topologies2 = determine_topology(mof14, structure=StructureType.MOF, clusterings=[Clustering.Auto, Clustering.Standard, Clustering.PE])
2 interpenetrated substructures:
⋅ Subnet 1 → AllNodes,SingleNodes,Standard: pto | PE: sqc11259
⋅ Subnet 2 → AllNodes,SingleNodes,Standard: pto | PE: sqc11259
julia> CrystalNets.total_interpenetration(topologies2, Clustering.AllNodes)
Dict{TopologicalGenome, Int64} with 1 entry:
pto => 2
julia> CrystalNets.total_interpenetration(topologies2)
Dict{TopologicalGenome, Int64} with 2 entries:
pto => 2
sqc11259 => 2