scona package¶
Module contents¶
scona¶
scona is a Python package for the analysis of structural covariance brain networks.
- Website (including documentation)::
- Source::
- Bug reports::
Simple example¶
FILL
Bugs¶
Please report any bugs that you find in an issue on our GitHub repository. Or, even better, fork the repository and create a pull request.
Submodules¶
scona.make_corr_matrices module¶
Tools to create a correlation matrix from regional measures
-
scona.make_corr_matrices.corrmat_from_regionalmeasures(regional_measures, names, covars=None, method='pearson')[source]¶ Calculate the correlation of names columns over the rows of regional_measures after correcting for covariance with the columns in covars
- Parameters
regional_measures (
pandas.DataFrame) – a pandas data frame with subjects as rows, and columns including brain regions and covariates. Should be numeric for the columns in names and covars_listnames (
list) – a list of the brain regions you wish to correlatecovars (
list) – covars is a list of covariates (as df column headings) to correct for before correlating the regions.methods (
string) – the method of correlation passed topandas.DataFrame.corr()
- Returns
A correlation matrix with rows and columns keyed by names
- Return type
-
scona.make_corr_matrices.create_corrmat(df_res, names=None, method='pearson')[source]¶ Correlate over the rows of df_res
- Parameters
df_res (
pandas.DataFrame) – df_res contains structural data about regions of the brain with subjects as rows after correction for any covariates of no interest.names (
list, optional) – The brain regions you wish to correlate over. These will become nodes in your graph. If names is None then all columns are included. Default is None.methods (
string, optional) – The method of correlation passed topandas.DataFrame.corr(). Default is pearsons correlation (pearson).
- Returns
A correlation matrix.
- Return type
- Raises
TypeError – If there are non numeric entries in the columns in df_res specified by names.
-
scona.make_corr_matrices.create_residuals_df(df, names, covars=[])[source]¶ Calculate residuals of columns specified by names, correcting for the columns in covars_list.
- Parameters
df (
pandas.DataFrame) – A pandas data frame with subjects as rows and columns including brain regions and covariates. Should be numeric for the columns in names and covars.names (
list) – A list of the brain regions you wish to correlate.covars (
list, optional) – A list of covariates to correct for before correlating the regional measures. Each element should correspond to a column heading in df. Default is an empty list.
- Returns
Residuals of columns names of df, correcting for covars
- Return type
- Raises
TypeError – if there are non numeric entries in the columns specified by names or covars
-
scona.make_corr_matrices.get_non_numeric_cols(df)[source]¶ returns the columns of df whose dtype is not numeric (e.g not a subtype of
numpy.number)- Parameters
df (
pandas.DataFrame) –- Returns
non numeric columns of df
- Return type
scona.make_graphs module¶
-
scona.make_graphs.anatomical_copy(G, nodal_keys=['name', 'name_34', 'name_68', 'hemi', 'centroids', 'x', 'y', 'z'], graph_keys=['parcellation', 'centroids'])[source]¶ Create a new graph from
Gpreserving: * nodes * edges * any nodal attributes specified in nodal_keys * any graph attributes specified in graph_keys- Parameters
G (
networkx.Graph) –Ghas anatomical data to copy toRnodal_keys (
list, optional) – The list of keys to treat as anatomical nodal data. Default is set toanatomical_node_attributes(), e.g["name", "name_34", "name_68", "hemi", "centroids", "x", "y", "z"]graph_keys (
list, optional) – The list of keys rto treat as anatomical graph data. Set toanatomical_graph_attributes()by default. E.g["centroids", "parcellation"]
- Returns
A new graph with the same nodes and edges as
Gand identical anatomical data.- Return type
See also
BrainNetwork.anatomical_copy(),copy_anatomical_data()
-
scona.make_graphs.anatomical_graph_attributes()[source]¶ default anatomical graph attributes for scona
- Returns
graph attributes considered “anatomical” by scona
- Return type
-
scona.make_graphs.anatomical_node_attributes()[source]¶ default anatomical nodal attributes for scona
- Returns
nodal attributes considered “anatomical” by scona
- Return type
-
scona.make_graphs.assign_node_centroids(G, centroids)[source]¶ Modify the nodal attributes “centroids”, “x”, “y”, and “z” of nodes of G, inplace. “centroids” will be set according to the scheme
G[i]["centroids"] = centroids[i]for nodesiinG. “x”, “y” and “z” are assigned as the first, second and third coordinate ofcentroids[i]respectively.- Parameters
G (
networkx.Graph) –centroids (
list) –centroids[i]is a tuple representing the cartesian coordinates of nodeiinG
- Returns
graph with nodal attributes modified
- Return type
-
scona.make_graphs.assign_node_names(G, parcellation)[source]¶ Modify nodal attribute “name” for nodes of G, inplace.
- Parameters
G (
networkx.Graph) –parcellation (
list) –parcellation[i]is the name of nodeiinG
- Returns
graph with nodal attributes modified
- Return type
-
scona.make_graphs.copy_anatomical_data(R, G, nodal_keys=['name', 'name_34', 'name_68', 'hemi', 'centroids', 'x', 'y', 'z'], graph_keys=['parcellation', 'centroids'])[source]¶ Copies nodal and graph attributes data from
GtoRfor keys included innodal_keysandgraph_keysrespectively.RandGare both graphs. If they have non equal vertex sets node data will only be copied fromGtoRfor nodes in the intersection.- Parameters
R (
networkx.Graph) –G (
networkx.Graph) –Ghas anatomical data to copy toRnodal_keys (
list, optional) – The list of keys to treat as anatomical nodal data. Default is set toanatomical_node_attributes(), e.g["name", "name_34", "name_68", "hemi", "centroids", "x", "y", "z"]graph_keys (
list, optional) – The list of keys rto treat as anatomical graph data. Set toanatomical_graph_attributes()by default. E.g["centroids", "parcellation"]
- Returns
graph
Rwith the anatomical data of graphG- Return type
See also
-
scona.make_graphs.get_random_graphs(G, n=10, Q=10, seed=None)[source]¶ Create n random graphs through edgeswapping.
- Parameters
G (
networkx.Graph) –n (
int, optional) –Q (
int, optional) – constant to specify how many swaps to conduct for each edge in Gseed (
int,random_stateorNone (default)) – Indicator of random state to pass to networkx
- Returns
A list of length n of randomisations of G created using
scona.make_graphs.random_graph()- Return type
list of
networkx.Graph
-
scona.make_graphs.graph_at_cost(M, cost, mst=True)[source]¶ Create a binary graph by thresholding weighted matrix M.
First creates the minimum spanning tree for the graph, and then adds in edges according to their connection strength up to cost.
- Parameters
M (
numpy.arrayorpandas.DataFrame) – A correlation matrix.cost (
float) – A number between 0 and 100. The resulting graph will have thecost*n/100highest weighted edges available, wherenis the number of edges in Gmst (
bool, optional) – IfFalse, skip creation of minimum spanning tree. This may cause output graph to be disconnected
- Returns
A binary graph
- Return type
- Raises
Exception – if M is not a
numpy.arrayorpandas.DataFrame
-
scona.make_graphs.is_anatomical_match(G, H, nodal_keys=['name', 'name_34', 'name_68', 'hemi', 'centroids', 'x', 'y', 'z'], graph_keys=['parcellation', 'centroids'])[source]¶ Check that G and H have identical anatomical data (including vertex sets).
- Parameters
G (
networkx.Graph) –H (
networkx.Graph) –nodal_keys (
list, optional) – The list of keys to treat as anatomical nodal data. Default is set toanatomical_node_attributes(), e.g["name", "name_34", "name_68", "hemi", "centroids", "x", "y", "z"]graph_keys (
list, optional) – The list of keys to treat as anatomical graph data. Set toanatomical_graph_attributes()by default. E.g["centroids", "parcellation"]
- Returns
Trueif G and H have the same anatomical data;Falseotherwise- Return type
-
scona.make_graphs.is_nodal_match(G, H, keys=None)[source]¶ Check that G and H have equal vertex sets. If keys is passed, also check that the nodal dictionaries of G and H (e.g the nodal attributes) are equal over the attributes in keys.
- Parameters
G (
networkx.Graph) –H (
networkx.Graph) –keys (
list, optional) – a list of attributes on which the nodal dictionaries of G and H should agree.
- Returns
Trueif G and H have equal vertex sets, or if keys is specified:Trueif vertex sets are equal AND the graphs’ nodal dictionaries agree on all attributes in keys.Falseotherwise- Return type
-
scona.make_graphs.random_graph(G, Q=10, seed=None)[source]¶ Return a connected random graph that preserves degree distribution by swapping pairs of edges, using
networkx.double_edge_swap().- Parameters
G (
networkx.Graph) –Q (
int, optional) – constant that specifies how many swaps to conduct for each edge in Gseed (
int,random_stateorNone (default)) – Indicator of random state to pass to networkx
- Returns
- Return type
- Raises
Exception – if it is not possible in 15 attempts to create a connected random graph.
-
scona.make_graphs.scale_weights(G, scalar=- 1, name='weight')[source]¶ Multiply edge weights of G by scalar.
- Parameters
G (
networkx.Graph) –scalar (
float, optional) – scalar value to multiply edge weights by. Default is -1name (
str, optional) – string that indexes edge weights. Default is “weight”
- Returns
- Return type
-
scona.make_graphs.threshold_graph(G, cost, mst=True)[source]¶ Create a binary graph by thresholding weighted graph G.
First creates the minimum spanning tree for the graph, and then adds in edges according to their connection strength up to cost.
- Parameters
G (
networkx.Graph) – A complete weighted graphcost (
float) – A number between 0 and 100. The resulting graph will have thecost*n/100highest weighted edges from G, wherenis the number of edges in G.mst (
bool, optional) – IfFalse, skip creation of minimum spanning tree. This may cause output graph to be disconnected
- Returns
A binary graph
- Return type
- Raises
Exception – If it is impossible to create a minimum_spanning_tree at the given cost
See also
scona.BrainNetwork.threshold()
-
scona.make_graphs.weighted_graph_from_df(df, create_using=None)[source]¶ Create a weighted graph from a correlation matrix.
- Parameters
df (
pandas.DataFrame) – a correlation matrixcreate_using (
networkx.Graph) – Use specified graph for result. The default is Graph()
- Returns
A weighted graph with edge weights equivalent to DataFrame entries
- Return type
-
scona.make_graphs.weighted_graph_from_matrix(M, create_using=None)[source]¶ Create a weighted graph from an correlation matrix.
- Parameters
M (
numpy.array) – an correlation matrixcreate_using (
networkx.Graph, optional) – Use specified graph for result. The default is Graph()
- Returns
A weighted graph with edge weights equivalent to matrix entries
- Return type
scona.graph_measures module¶
-
scona.graph_measures.assign_interhem(G)[source]¶ Assigns nodal and edge attributes of G. Modifies G in place. An edge is considered interhemispheric if the x coordinates of its nodes have different signs.
Edge attributes:
- “interhem”int
1 if the edge is interhemispheric, 0 otherwise
Node attributes
- “hemisphere”str
L or R, as determined by the sign of the x coordinate and assuming MNI space. The x coordinates are negative in the left hemisphere and positive in the right.
- “interhem”int
the number of adjacent interhemispheric edges
- “interhem_proportion”float
the proportion of adjacent edges that are interhemispheric
- Parameters
G (
networkx.Graph) – a graph with nodal attribute ‘centroids’ or ‘x’ defined for each node. The value of ‘centroids’ should be the cartesian coordinates of each node.- Returns
G
- Return type
-
scona.graph_measures.assign_nodal_distance(G)[source]¶ Assigns nodal and edge attributes of G. Modifies G in place.
Edge attributes
- “euclidean”float
the euclidean length, derived from the “centroids” values of nodes
Node attributes
- “total_dist”float
the total length of the incident edges
- “average_dist”float
the average length of the incident edges
- Parameters
G (
networkx.Graph) – a graph with nodal attribute ‘centroids’ defined for each node. The value of ‘centroids’ should be the cartesian coordinates of each node.- Returns
G
- Return type
-
scona.graph_measures.calc_modularity(G, nodal_partition)[source]¶ Calculate the modularity of G under partition nodal_partition.
- Parameters
G (
networkx.Graph) –nodal_partition (
dict) – a dictionary nodes to communities
- Returns
the modularity of G
- Return type
-
scona.graph_measures.calc_nodal_partition(G)[source]¶ Calculate a nodal partition of G using the louvain algorithm as iBrainNetworkommunity.best_partition`
Note that this is a time intensive process and it is also non-deterministic, so for consistency and speed it’s best to hold on to your partition.
- Parameters
G (
networkx.Graph) – A binary graph- Returns
Two dictionaries represent the resulting nodal partition of G. The first maps nodes to modules and the second maps modules to nodes.
- Return type
-
scona.graph_measures.calculate_global_measures(G, partition=None, existing_global_measures=None)[source]¶ Calculate global measures average_clustering, average_shortest_path_length, assortativity, modularity, and efficiency of G.
Note: Global measures will not be calculated again if they have already been calculated. So it is only needed to calculate them once and then they aren’t calculated again.
- Parameters
G (
networkx.Graph) – A binary graphpartition (
dict, optional) – A nodal partition of G. A dictionary mapping nodes of G to modules. Pass a partition in order to calculate the modularity of G.existing_global_measures (
dict, optional) – An existing dictionary of global measures of G can be passed.calculate_global_measures()will not recalculate any measures already indexed in G
- Returns
a dictionary of global network measures of G
- Return type
See also
scona.BrainNetwork.calculate_global_measures()
-
scona.graph_measures.calculate_nodal_measures(G, partition=None, measure_list=None, additional_measures=None, force=True)[source]¶ Calculate and store nodal measures as nodal attributes.
By default calculate_nodal_measures calculates the following :
“degree” : int
“closeness” : float
“betweenness” : float
“shortest_path_length” : float
“clustering” : float
“participation_coefficient” : float
Use measure_list to specify which of the default nodal attributes to calculate. Use additional_measures to describe and calculate new measure definitions.
- Parameters
G (
networkx.Graph) –measure_list (
listofstr, optional) – pass a subset of of the keys defined above to specify which of the default measures to calculateadditional_measures (
dict, optional) – map from names of nodal attributes to functions defining how they should be calculated. Such a function should take a graph as an argument and return a dictionary mapping nodes to attribute values.force (
bool, optional) – pass True to recalculate any measures that already exist in the nodal attributes.
See also
BrainNetwork.calculate_nodal_measures(),calc_nodal_partition()Example
-
scona.graph_measures.participation_coefficient(G, module_partition)[source]¶ Computes the participation coefficient of nodes of G with partition defined by module_partition. (Guimera et al. 2005).
- Parameters
G (
networkx.Graph) –module_partition (
dict) – a dictionary mapping each community name to a list of nodes in G
- Returns
a dictionary mapping the nodes of G to their participation coefficient under the participation specified by module_partition.
- Return type
-
scona.graph_measures.rich_club(G)[source]¶ Calculate the rich club coefficient of G for each degree between 0 and
max([degree(v) for v in G.nodes]).- Parameters
G (
networkx.Graph) – a binary graph- Returns
a dictionary mapping integer
xto the rich club coefficient of G for degreex- Return type
See also
BrainNetwork.rich_club()
-
scona.graph_measures.shortest_path(G)[source]¶ Calculate average shortest path length for each node in G. “length” in this case means the number of edges, and does not consider euclidean distance.
- Parameters
G (
networkx.Graph) – a connected graph- Returns
a dictionary mapping a node v to the average length of the shortest from v to other nodes in G.
- Return type
-
scona.graph_measures.small_world_coefficient(G, R)[source]¶ Calculate the small world coefficient of G relative to R.
Small coefficient is (G.average_clustering/R.average_clustering) / (G.average_shortest_path_length / R.average_shortest_path_length) , where average_clustering and average_shortest_path_length are a graph’s global measures.
- Parameters
G (
networkx.Graph) – A binary graphR (
networkx.Graph) – A binary graph
- Returns
The small world coefficient of G relative to R
- Return type
-
scona.graph_measures.small_world_sigma(tupleG, tupleR)[source]¶ Compute small world sigma from tuples
-
scona.graph_measures.z_score(G, module_partition)[source]¶ Calculate the z-score of the nodes of G under partition module_partition.
- Parameters
G (
networkx.Graph) –module_partition (
dict) – a dictionary mapping each community name to a lists of nodes in G
- Returns
a dictionary mapping the nodes of G to their z-score under module_partition.
- Return type
scona.classes module¶
-
class
scona.classes.BrainNetwork(network=None, parcellation=None, centroids=None)[source]¶ Bases:
networkx.classes.graph.GraphA lightweight subclass of
networkx.Graph.- Parameters
network (
networkx.Graphornumpy.ndarrayorpandas.DataFrame, optional) – network is used to define graph features of G. If an array is passed, create a weighted graph from this array. If a graph is passed, G will be aBrainNetworkotherwise identical to this graph.parcellation (
listofstr, optional) – A list of node names, passed toBrainNetwork.set_parcellation()centroids (
listoftuple, optional) – A list of node centroids, passed toBrainNetwork.set_centroids()
-
anatomical_graph_attributes¶ a list of graph attribute keys to treat as anatomical.
See also
-
anatomical_copy()[source]¶ Create a new graph from G preserving: * nodes * edges * any nodal attributes specified in G.anatomical_node_attributes * any graph attributes specified in G.anatomical_graph_attributes *
G.anatomical_node_attributes*G.anatomical_graph_attributes- Returns
A new graph with the same nodes and edges as G and identical anatomical data.
- Return type
See also
BrainNetwork.anatomical_copy(),anatomical_data(),set_anatomical_node_attributes(),set_anatomical_graph_attributes()
-
calculate_global_measures(force=False, seed=None, partition=True)[source]¶ Calculate global measures average_clustering, average_shortest_path_length, assortativity, modularity, and efficiency of G. The resulting dictionary of global measures can be accessed and manipulated as
G.graph['global_measures']- Parameters
force (
bool) – pass True to recalculate any global measures that have already been calculatedpartition (
bool) – The “modularity” measure evaluates a graph partition. pass True to calculate the partition of each graph usingBrainNetwork.partition(). Note that this won’t recalculate a partition that exists. If False, modularity may not be calculated.
- Returns
a dictionary of global network measures of G
- Return type
See also
-
calculate_nodal_measures(force=False, measure_list=None, additional_measures=None)[source]¶ Calculate and store nodal measures as nodal attributes which can be accessed directly, or using
BrainNetwork.report_nodal_measures()`()By default calculates:
betweenness : float
closeness : float
clustering : float
degree : int
module : int
participation_coefficient : float
shortest_path_length : float
Use measure_list to specify which of the default nodal attributes to calculate. Use additional_measures to describe and calculate new measure definitions.
- Parameters
measure_list (
listofstr, optional) – pass a subset of of the keys defined above to specify which of the default measures to calculateadditional_measures (
dict, optional) – map from names of nodal attributes to functions defining how they should be calculated. Such a function should take a graph as an argument and return a dictionary mapping nodes to attribute values.force (
bool, optional) – pass True to recalculate any measures that already exist in the nodal attributes.
Example
-
calculate_spatial_measures()[source]¶ Assigns spatial nodal and edge attributes of G. Modifies G in place.
Edge attributes:
- “euclidean”float
the euclidean length, as derived from the nodal attribute “centroids”
- “interhem”int
1 if the edge is interhemispheric, 0 otherwise
Node attributes
- “total_dist”float
the total length of the incident edges
- “average_dist”float
the average length of the incident edges
- “hemisphere”str
L or R, as determined by the sign of the x coordinate and assuming MNI space. The x coordinates are negative in the left hemisphere and positive in the right.
- “interhem”int
the number of adjacent interhemispheric edges
- “interhem_proportion”float
the proportion of adjacent edges that are interhemispheric
- Raises
KeyError – If the graph centroids are not set
See also
BrainNetwork.set_centroids(),assign_interhem(),assign_nodal_distance()
-
partition(force=False)[source]¶ Calculate the best nodal partition of G using the louvain algorithm as implemented in
community.best_partition(). If desired, the node-wise partition can be accessed and manipulated by the nodal attribute “module” and the module-wise partition asG.graph["partition"]- Parameters
force (
bool) – pass True to recalculate when a partition exists- Returns
Two dictionaries represent the resulting nodal partition of G. The first maps nodes to modules and the second maps modules to nodes.
- Return type
See also
BrainNetwork.calc_nodal_partition()
-
report_nodal_measures(columns=None, as_dict=False)[source]¶ Report the nodal attributes of G as a pandas dataframe or python dictionary
- Parameters
- Returns
the node attribute data from G as a pandas dataframe or dict of dicts
- Return type
pandas.DataFrameor dict
-
rich_club(force=False)[source]¶ Calculate the rich club coefficient of G for each degree between 0 and
max([degree(v) for v in G.nodes]). The resulting dictionary of rich club coefficients can be accessed and manipulated asG.graph['rich_club']- Parameters
force (
bool) – pass True to recalculate when a dictionary of rich club coefficients already exists- Returns
a dictionary mapping integer x to the rich club coefficient of G for degree x
- Return type
See also
-
set_anatomical_graph_attributes(names)[source]¶ Define the list of graph attribute keys to preserve when using
BrainNetwork.anatomical_copy(). This list is set usinganatomical_graph_attributes()when a BrainNetwork is initialised, and can be accessed asG.anatomical_graph_attributesIt is also among the object attributes preserved byBrainNetwork.anatomical_copy()- Parameters
names (
list) – a list of graph attribute keys to treat as anatomical.
-
set_anatomical_node_attributes(names)[source]¶ Define the list of node attribute keys to preserve when using
BrainNetwork.anatomical_copy(). This list is set usinganatomical_node_attributes()when a BrainNetwork is initialised, and can be accessed asG.anatomical_node_attributesIt is also among the object attributes preserved byBrainNetwork.anatomical_copy()- Parameters
names (
list) – a list of node attribute keys to treat as anatomical.
-
set_centroids(centroids)[source]¶ Modify the nodal attributes “centroids”, “x”, “y”, and “z” of nodes of G inplace.
- Parameters
G (
networkx.Graph) –centroids (
list) –centroids[i]is a tuple representing the cartesian coordinates of nodeiinG. “x”, “y” and “z” will be assigned as the first, second and third coordinate ofcentroids[i]respectively.
See also
BrainNetwork.calculate_spatial_measures(),BrainNetwork.set_parcellation(),assign_node_centroids()
-
set_parcellation(parcellation)[source]¶ Modify nodal attribute “name” for nodes of G inplace.
- Parameters
G (
networkx.Graph) –parcellation (
list) –parcellation[i]is the name of nodeiinG
See also
assign_node_names(),BrainNetwork.set_centroids()
-
threshold(cost, mst=True)[source]¶ Create a binary graph by thresholding weighted BrainNetwork G
First creates the minimum spanning tree for the graph, and then adds in edges according to their connection strength up to cost.
- Parameters
- Returns
A binary graph
- Return type
- Raises
Exception – If it is impossible to create a minimum_spanning_tree at the given cost
See also
threshold_graph()
-
class
scona.classes.GraphBundle(graph_list, name_list)[source]¶ Bases:
dictThe GraphBundle class (after instantiating - object) is the scona way to handle across-network comparisons. What is it? Essentially it’s a python dictionary with BrainNetwork objects as values (
str:BrainNetworkpairs).Mainly used to create random graphs for comparison with your original network data.
- Parameters
graph_list (
listofnetworkx.Graph) –
See also
Example
-
add_graphs(graph_list, name_list=None)[source]¶ Update dictionary with graph_list : names_list pairs.
- Parameters
graph_list (
listofnetworkx.Graph) –
See also
-
anatomical_matches(nodal_keys=['name', 'name_34', 'name_68', 'hemi', 'centroids', 'x', 'y', 'z'], graph_keys=['parcellation', 'centroids'])[source]¶ Checks that all graphs in GraphBundle are pairwise anatomical matches as defined in
is_anatomical_match().- Parameters
- Returns
True if all graphs are anatomically matched, False otherwise.
- Return type
See also
is_anatomical_match(),BrainNetwork.is_nodal_match()
-
apply(graph_function)[source]¶ Apply a user defined function to each network in a
GraphBundle.- Parameters
graph_function (
types.FunctionType) – Function defined on aBrainNetworkobject
Example
To calculate and return the degree for each network in a
GraphBundlepass this following expression to apply as the graph_function.
-
create_random_graphs(gname, n, Q=10, name_list=None, rname='_R', seed=None)[source]¶ Create n edge swap randomisations of
BrainNetworkkeyed by gname. These random graphs are added to GraphBundle.- Parameters
gname (
str) – indexes a graph in GraphBundlen (
int) – the number of random graphs to createQ (
int, optional) – constant to specify how many swaps to conduct for each edge in Gname_list (
listofstr, optional) – a list of names to use for indexing the new random graphs in GraphBundle.rname (
str, optional) – ifname_list=Nonethe new random graphs will be indexed according to the schemegname + rname + rwhere r is some integer.seed (
int,random_stateorNone (default)) – Indicator of random state to pass tonetworkx.double_edge_swap()
See also
get_random_graphs(),random_graph(),BrainNetwork.add_graphs()
-
nodal_matches()[source]¶ Checks the statement “All graphs in GraphBundle have congruent vertex sets”
- Returns
True if all graphs have the same node set, False otherwise.
- Return type
See also
is_nodal_match(),BrainNetwork.anatomical_matches()
-
report_global_measures(as_dict=False, partition=True)[source]¶ Calculate global_measures for each BrainNetwork object and report as a
pandas.DataFrameor nested dict.Note: Global measures will not be calculated again if they have already been calculated. So it is only needed to calculate them once and then they aren’t calculated again.
- Parameters
as_dict (
bool) – pass True to return global measures as a nested dictionary; pass False to return apandas.DataFramepartition (
bool) – argument to pass toBrainNetwork.calculate_global_measures()
- Returns
- Return type
pandas.DataFrameor dict
-
report_rich_club(as_dict=False)[source]¶ Calculate rich_club coefficients for each BrainNetwork object and report as a
pandas.DataFrameor nested dict.- Parameters
as_dict (
bool) – pass True to return rich club coefficients as a nested dictionary; pass False to return apandas.DataFrame- Returns
- Return type
pandas.DataFrameor dict
See also
-
report_small_world(gname)[source]¶ Calculate the small world coefficient of gname relative to each other graph in GraphBundle.
- Parameters
gname (
str) – indexes a graph in GraphBundle- Returns
a dictionary, mapping a GraphBundle key “x” to the small coefficient of graph “gname” relative to graph “x”.
- Return type
See also
small_world_coefficient()
scona.stats_functions module¶
-
scona.stats_functions.residuals(x, y)[source]¶ Return residuals of least squares solution to y = AB where A =[[x 1]]. Uses
numpy.linalg.lstsq()to find B- Parameters
x (
array) –y (
array) –
- Returns
residuals of least squares solution to y = [[x 1]]B
- Return type
array