Classes and Functions

Graph classes

class NEMtropy.graph_classes.DirectedGraph(adjacency=None, edgelist=None, degree_sequence=None, strength_sequence=None)[source]

Directed graph instance can be initialised with adjacency matrix, edgelist, degree sequence or strengths sequence.

Parameters
  • adjacency (numpy.ndarray, list, scipy.sparse_matrix, optional) – Adjacency matrix, defaults to None.

  • edgelist (numpy.ndarray, list, optional) – edgelist, defaults to None.

  • degree_sequence (numpy.ndarray, optional) – degrees sequence, defaults to None.

  • strength_sequence (numpy.ndarray, optional) – strengths sequence, defaults to None.

clean_edges()[source]

Deletes all initialized graph attributes

degree_reduction()[source]

Carries out degree reduction for DBCM. The graph should be initialized.

ensemble_sampler(n, cpu_n=1, output_dir='sample/', seed=42)[source]

The function sample a given number of graphs in the ensemble generated from the last model solved. Each grpah is an edgelist written in the output directory as .txt file. The function is parallelised and can run on multiple cpus.

Parameters
  • n (int) – Number of graphs to sample.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 1.

  • output_dir (str, optional) – Name of the output directory, defaults to “sample/”.

  • seed (int, optional) – Random seed, defaults to 42.

Raises

ValueError – [description]

set_adjacency_matrix(adjacency)[source]

Initializes a graph from the adjacency matrix.

Parameters

adjacency (numpy.ndarray, list, scipy.sparse_matrix) – Adjacency matrix.

set_degree_sequences(degree_sequence)[source]

Initializes graph from the degrees sequence.

Parameters

adjacency (numpy.ndarray) – Degrees sequence

set_edgelist(edgelist)[source]

Initializes a graph from the edgelist.

Parameters

adjacency (numpy.ndarray, list) – Edgelist

solve_tool(model, method, initial_guess=None, adjacency=None, method_adjacency='newton', initial_guess_adjacency='random', max_steps=100, full_return=False, verbose=False, linsearch=True, tol=1e-08, eps=1e-08)[source]

The function solves the ERGM optimization problem from a range of available models. The user can choose among three optimization methods. The graph should be initialized.

Parameters
  • model (str) – Available models are: - dcm: solves DBCM respect to the parameters x and “y” of the loglikelihood function, it works for uweighted directed graphs [insert ref]. - dcm_exp: differently from the dcm option, dcm_exp considers the exponents of x and y as parameters [insert ref]. - decm: solves DECM respect to the parameters a_out, a_in, b_out and b_in of the loglikelihood function, it is conceived for weighted directed graphs [insert ref]. - decm_exp: differently from the decm option, decm_exp considers the exponents of a_out, a_in, b_out and b_in* as parameters [insert ref]. - crema: solves CReMa for a weighted directd graphs. In order to compute beta parameters, it requires information about the binary structure of the network. These can be provided by the user by using adjacency paramenter. - crema-sparse: alternative implementetio of crema for large graphs. The creama-sparse model doesn’t compute the binary probability matrix avoing memory problems for large graphs.

  • method (str) – Available methods to solve the given model are: - newton: uses Newton-Rhapson method to solve the selected model, it can be memory demanding for crema because it requires the computation of the entire Hessian matrix. This method is not available for creama-sparse. - quasinewton: uses Newton-Rhapson method with Hessian matrix approximated by its principal diagonal to find parameters maximising loglikelihood function. - fixed-point: uses a fixed-point method to find parameters maximising loglikelihood function.

  • initial_guess (str, optional) –

    Starting point solution may affect the results of the optization process. The user can provid an initial guess or choose between the following options:

    • Binary Models:

      • random: random numbers in (0, 1);

      • uniform: uniform initial guess in (0, 1);

      • degrees: initial guess of each node is proportianal to its degree;

      • degrees_minor: initial guess of each node is inversely proportional to its degree;

      • chung_lu: initial guess given by Chung-Lu formula;

    • Weighted Models:

      • random: random numbers in (0, 1);

      • uniform: uniform initial guess in (0, 1);

      • strengths: initial guess of each node is proportianal to its stength;

      • strengths_minor: initial guess of each node is inversely proportional to its strength;

  • adjacency (str or numpy.ndarray, optional) – Adjacency can be a binary method (defaults is dcm_exp) or an adjacency matrix.

  • method_adjacency (str, optional) – If adjacency is a model, it is the methdod used to solve it. Defaults to “newton”.

  • initial_guess_adjacency (str, optional) – If adjacency is a model, it is the chosen initial guess. Defaults to “random”.

  • max_steps (int, optional) – maximum number of iteration, defaults to 100.

  • full_return (bool, optional) – If True the algorithm returns more statistics than the obtained solution, defaults to False.

  • verbose (bool, optional) – If True the algorithm prints a bunch of statistics at each step, defaults to False.

  • linsearch (bool, optional) – If True the linsearch function is active, defaults to True.

  • tol (float, optional) – parameter controlling the tollerance of the norm the gradient function, defaults to 1e-8.

  • eps (float, optional) – parameter controlling the tollerance of the difference between two iterations, defaults to 1e-8.

class NEMtropy.graph_classes.UndirectedGraph(adjacency=None, edgelist=None, degree_sequence=None, strength_sequence=None)[source]

Undirected Graph instance can be initialised with adjacency matrix, edgelist, degrees sequence or strengths sequence.

Parameters
  • adjacency (numpy.ndarray, list or scipy.sparse_matrix, optional) – Adjacency matrix, defaults to None.

  • edgelist (numpy.ndarray, list, optional) – edgelist, defaults to None.

  • degree_sequence (numpy.ndarray, optional) – degrees sequence, defaults to None.

  • strength_sequence (numpy.ndarray, optional) – strengths sequence, defaults to None.

clean_edges()[source]

Deletes all the initialiased attributes.

degree_reduction()[source]

Carries out degree reduction for UBCM. The graph should be initialized.

ensemble_sampler(n, cpu_n=1, output_dir='sample/', seed=42)[source]

The function sample a given number of graphs in the ensemble generated from the last model solved. Each grpah is an edgelist written in the output directory as .txt file. The function is parallelised and can run on multiple cpus.

Parameters
  • n (int) – Number of graphs to sample.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 1.

  • output_dir (str, optional) – Name of the output directory, defaults to “sample/”.

  • seed (int, optional) – Random seed, defaults to 42.

Raises

ValueError – [description]

set_adjacency_matrix(adjacency)[source]

Initialises graph given the adjacency matrix.

Parameters

adjacency (numpy.ndarray, list, scipy.sparse_matrix) – ajdacency matrix.

set_degree_sequences(degree_sequence)[source]

Initialises graph given the degrees sequence.

Parameters

degree_sequence (numpy.ndarray) – degrees sequence.

set_edgelist(edgelist)[source]

Initialises graph given the edgelist.

Parameters

edgelist (numpy.ndarray, list) – edgelist.

solve_tool(model, method, initial_guess=None, adjacency='cm_exp', method_adjacency='newton', initial_guess_adjacency='random', max_steps=100, full_return=False, verbose=False, linsearch=True, tol=1e-08, eps=1e-08)[source]

The function solves the ERGM optimization problem from a range of available models. The user can choose among three optimization methods. The graph should be initialized.

Parameters
  • model (str) –

    Available models are:

    • cm: solves UBCM respect to the parameters x of the mof.loglikelihood function, it works for uweighted undirected graphs [insert ref].

    • cm_exp: differently from the cm option, cm_exp considers the exponents of x as parameters [insert ref].

    • ecm: solves UECM respect to the parameters x and y of the mof.loglikelihood function, it is conceived for weighted undirected graphs [insert ref].

    • ecm_exp: differently from the ecm option, ecm_exp considers the exponents of x and y as parameters [insert ref].

    • crema: solves CReMa for a weighted undirectd graphs. In order to compute beta parameters, it requires information about the binary structure of the network. These can be provided by the user by using adjacency paramenter.

    • crema-sparse: alternative implementetio of crema for large graphs. The creama-sparse model doesn’t compute the binary probability matrix avoing memory problems for large graphs.

  • method (str) –

    Available methods to solve the given model are:

    • newton: uses Newton-Rhapson method to solve the selected model, it can be memory demanding for crema because it requires the computation of the entire Hessian matrix. This method is not available for creama-sparse.

    • quasinewton: uses Newton-Rhapson method with Hessian matrix approximated by its principal diagonal to find parameters maximising mof.loglikelihood function.

    • fixed-point: uses a fixed-point method to find parameters maximising mof.loglikelihood function.

  • initial_guess (str, optional) –

    Starting point solution may affect the results of the optization process. The user can provid an initial guess or choose between the following options:

    • Binary Models:
      • random: random numbers in (0, 1);

      • uniform: uniform initial guess in (0, 1);

      • degrees: initial guess of each node is proportianal to its degree;

      • degrees_minor: initial guess of each node is inversely proportional to its degree;

      • chung_lu: initial guess given by Chung-Lu formula;

    • Weighted Models:
      • random: random numbers in (0, 1);

      • uniform: uniform initial guess in (0, 1);

      • strengths: initial guess of each node is proportianal to its stength;

      • strengths_minor: initial guess of each node is inversely proportional to its strength;

  • adjacency (str or numpy.ndarray, optional) – Adjacency can be a binary method (defaults is cm_exp) or an adjacency matrix.

  • method_adjacency (str, optional) – If adjacency is a model, it is the methdod used to solve it. Defaults to “newton”.

  • initial_guess_adjacency (str, optional) – If adjacency is a model, it is the chosen initial guess. Defaults to “random”.

  • max_steps (int, optional) – maximum number of iteration, defaults to 100.

  • full_return (bool, optional) – If True the algorithm returns more statistics than the obtained solution, defaults to False.

  • verbose (bool, optional) – If True the algorithm prints a bunch of statistics at each step, defaults to False.

  • linsearch (bool, optional) – If True the linsearch function is active, defaults to True.

  • tol (float, optional) – parameter controlling the tollerance of the norm the gradient function, defaults to 1e-8.

  • eps (float, optional) – parameter controlling the tollerance of the difference between two iterations, defaults to 1e-8.

Ensemble generator

NEMtropy.ensemble_generator.ensemble_sampler_cm_graph(outfile_name, x, cpu_n=2, seed=None)[source]
Produce a single undirected binary graph after UBCM.

The graph is saved as an edges list on a file.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • x (numpy.ndarray) – UBCM solution

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_decm_det_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single directed weighted graph after CReMa,

given a fixed, complete network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Edges list: out-nodes array, in-nodes array and links weights.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_decm_prob_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single directed weighted graph after CReMa,

given a probabilistic network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Edges list: out-nodes array, in-nodes array and links probabilities.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_ecm_det_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single undirected weighted graph after CReMa,

given a fixed, complete network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Edges list: out-nodes array, in-nodes array and links weights.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_ecm_prob_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single undirected weighted graph after CReMa,

given a probabilistic network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Edges list: out-nodes array, in-nodes array and links probabilities.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_sparse_decm_prob_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single directed weighted graph after the

sparse version of CReMa, given a probabilistic network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray)) – Tuple of DBCM out-degrees parameters and in-degrees parameters solution.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_crema_sparse_ecm_prob_graph(outfile_name, beta, adj, cpu_n=2, seed=None)[source]
Produce a single undirected weighted graph after CReMa,

given a probabilistic network topology. The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • beta (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Edges list: out-nodes array, in-nodes array and links probabilities.

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_dcm_graph(outfile_name, x, y, cpu_n=2, seed=None)[source]
Produce a single directed binary graph after DBCM.

The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • x (numpy.ndarray) – DBCM out-degrees parameters solution

  • y (numpy.ndarray) – DBCM in-degrees parameters solution

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_decm_graph(outfile_name, a_out, a_in, b_out, b_in, cpu_n=2, seed=None)[source]
Produce a single directed weighted graph after DECM.

The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • a_out (numpy.ndarray) – DECM out-degrees parameters solution

  • a_in (numpy.ndarray) – DECM in-degrees parameters solution

  • b_out (numpy.ndarray) – DECM out-strengths parameters solution

  • b_in (numpy.ndarray) – DECM in-strengths parameters solution

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

NEMtropy.ensemble_generator.ensemble_sampler_ecm_graph(outfile_name, x, y, cpu_n=2, seed=None)[source]
Produce a single undirected weighted graph after ECM.

The graph is saved as an edges list on a file. The function is parallelyzed.

Parameters
  • outfile_name (str) – Name of the file where to save the graph on.

  • x (numpy.ndarray) – ECM degrees solution

  • y (numpy.ndarray) – ECM strengths solution

  • cpu_n (int, optional) – Number of cpus to use, defaults to 2

  • seed (int, optional) – Random seed, defaults to None

Returns

The outfile_name input parameter

Return type

str

The function randomly returns an undirected link

between two given nodes after the UBCM.

Parameters
  • args_1 ((int, float)) – Tuple containing node 1 and related parameter solution.

  • args_2 ((int, float)) – Tuple containing node 2 and related parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the couple of nodes, otherwise returns None

Return type

(int, int)

The function randomly returns a directed link with related weight
between two given nodes after the CReMa undirected for a deterministic

topology.

Parameters
  • args_1 ((int, float)) – Tuple containing out node and out-strength parameter solution.

  • args_2 ((int, float)) – Tuple containing in node and in-strength parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns a directed link with related weight

between two given nodes after the CReMa directed for a probabilistic topology.

Parameters
  • args_1 ((int, float)) – Tuple containing out node and out-strength parameter solution.

  • args_2 ((int, float)) – Tuple containing in node and in-strength parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns an undirected link with related weight

between two given nodes after the CReMa undirected for a deterministic topology.

Parameters
  • args_1 ((int, float)) – Tuple containing node 1, and strength parameter solution.

  • args_2 ((int, float)) – Tuple containing node 2, and strength parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns an undirected link with related weight

between two given nodes after the CReMa undirected for probabilistic o topology.

Parameters
  • args_1 ((int, float)) – Tuple containing node 1, and strength parameter solution.

  • args_2 ((int, float)) – Tuple containing node 2, and strength parameter solution.

  • p_ensemble (float) – Probability the link exists

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns a directed link with related weight

between two given nodes after the CReMa directed for a probabilistic topology. Function for the sparse version of CReMa.

Parameters
  • args_1 ((int, float, float)) – Tuple containing out node, out-strength parameter solution and out-degree parameter solution from the probabilistic model

  • args_2 ((int, float, float)) – Tuple containing in node and in-strength parameter solution and in-degree parameter solution from the probabilistic model

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns an undirected link with related weight

between two given nodes after the CReMa undirected for probabilistic o topology. Function for the sparse version of CReMa.

Parameters
  • args_1 ((int, float)) – Tuple containing node 1, strength parameter solution and degree parameter solution from the probabilistic model.

  • args_2 ((int, float)) – Tuple containing node 2, strength parameter solution and degree parameter solution from the probabilistic model

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns a directed link

between two given nodes after the DBCM.

Parameters
  • args_1 ((int, float)) – Tuple containing out node and out-degree parameter solution.

  • args_2 ((int, float)) – Tuple containing in node and in-degree parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the couple of nodes, otherwise returns None

Return type

(int, int)

The function randomly returns an directed link with related weight

between two given nodes after the DECM.

Parameters
  • args_1 ((int, float, float)) – Tuple containing out node, out-degree parameter solution and out-strength parameter solution.

  • args_2 ((int, float, float)) – Tuple containing in node, in-degree parameter solution and in-strength parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

The function randomly returns an undirected link with related weight

between two given nodes after the ECM.

Parameters
  • args_1 ((int, float, float)) – Tuple containing node 1, degree parameter solution and strength parameter solution.

  • args_2 ((int, float, float)) – Tuple containing node 2, degree parameter solution and strength parameter solution.

  • seed (int, optional) – Random seed, defaults to None

Returns

If the links exists returns the triplet of nodes and weight, otherwise returns None

Return type

(int, int, float)

Matrix generator

NEMtropy.matrix_generator.barabasi_albert_graph_nx(n, m, sup_ext, alpha, seed=None, is_weighted=None, is_sparse=False)[source]

Generates a undirected weighted adjacency matrix using barabasi-albert model for the binary part and assigning weights extracted from a uniform, gaussian or powerlaw distribution.

Parameters
  • n (int) – Number of nodes.

  • m (int) – Number of edges to attach from a new node to existing nodes.

  • sup_ext (int, float, optional) – Maximum strength, defaults to 10.

  • alpha (float) – Powerlaw exponent.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • is_weighted (bool, optional) – If True the adjacency matrix is weighted, defaults to None.

  • is_sparse (bool, optional) – If True the adjacency matrix is sparse, defaults to False.

Returns

Adjacency matrix.

Return type

numpy.ndarray, scipy.sparse_matrix

NEMtropy.matrix_generator.random_binary_matrix_generator_custom_density(n, p=0.1, sym=False, seed=None)[source]

Generates a random binary adjacency matrix with density determined by the parameter p.

Parameters
  • n (int) – Number of nodes.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

Returns

Binary adjacency matrix.

Return type

numpy.ndarray

NEMtropy.matrix_generator.random_binary_matrix_generator_custom_density_sparse(n, p=0.1, sym=False, seed=None)[source]

Generates a random binary adjacency sparse matrix with density determined by the parameter p.

Parameters
  • n (int) – Number of nodes.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – seed of the random process, defaults to None.

Returns

a binary sparse adjacency matrix.

Return type

scipy.sparse_matrix

NEMtropy.matrix_generator.random_binary_matrix_generator_dense(n, sym=False, seed=None)[source]

Generates a densely connected random binary adjacency matrix.

Parameters
  • n (int) – Number of nodes.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

Returns

Binary adjacency matrix.

Return type

numpy.ndarray

NEMtropy.matrix_generator.random_gaussian_weighted_matrix_generator_custom_density_sparse(n, mean, sigma, p=0.1, sym=False, seed=None, intweights=False)[source]

Generates a random weighted adjacency sparse matrix with density determined by the parameter p and weights extracted from a gaussian distribution.

Parameters
  • n (int) – Number of nodes.

  • mean (int, float) – Mean of the gaussian.

  • sigma (int, float) – Sigma of the gaussian.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sym (bool, optional) – If True the adjacency matrix is simmetrix, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • intweights (bool, optional) – If True weights are integers, defaults to False.

Returns

Weighted sparse adjacency matrix.

Return type

scipy.sparse_matrix

NEMtropy.matrix_generator.random_graph_nx(n, p, sup_ext, alpha, seed=None, is_weighted=None, is_sparse=False)[source]

Retunrs an undirected weighted adjacency matrix using erdos-renyi model to generate the binary adjacency and assigning weights extracted’ from a uniform, gaussian or powerlaw distribution.

Parameters
  • n (int) – umber of nodes.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sup_ext (int, float, optional) – Maximum strength, defaults to 10.

  • alpha (float) – Powerlaw exponent.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • is_weighted (bool, optional) – If True the adjacency matrix is weighted, defaults to None.

  • is_sparse (bool, optional) – If True the adjacency matrix is sparse, defaults to False.

Returns

Adjacency matrix.

Return type

numpy.ndarray, scipy.sparse_matrix

NEMtropy.matrix_generator.random_uniform_weighted_matrix_generator_custom_density_sparse(n, sup_ext, p=0.1, sym=False, seed=None, intweights=False)[source]

Generates a random weighted adjacency sparse matrix with density determined by the parameter p and weights extracted from a uniform distribution.

Parameters
  • n (int) – Number of nodes.

  • sup_ext (int, float, optional) – Maximum strength, defaults to 10.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sym (bool, optional) – If True the adjacency matrix is simmetrix, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • intweights (bool, optional) – If True weights are integers, defaults to False.

Returns

Weighted sparse adjacency matrix.

Return type

scipy.sparse_matrix

NEMtropy.matrix_generator.random_weighted_matrix_generator_dense(n, sup_ext=10, sym=False, seed=None, intweights=False)[source]

Generates a densely connected random weighted adjacency matrix.

Parameters
  • n (int) – Number of nodes.

  • sup_ext (int, float, optional) – Maximum strength, defaults to 10.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • intweights (bool, optional) – If True weights are integers, defaults to False.

Returns

Weighted adjacency matrix.

Return type

numpy.ndarray

NEMtropy.matrix_generator.random_weighted_matrix_generator_gaussian_custom_density(n, mean, sigma, p=0.1, sym=False, seed=None, intweights=False)[source]

Generates a random weighted adjacency matrix with density determined by the parameter p and weights extracted from a gaussian distribution.

Parameters
  • n (int) – Number of nodes.

  • mean (int, float) – Mean of the gaussian.

  • sigma (int, float) – Sigma of the gaussian.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • intweights (bool, optional) – If True weights are integers, defaults to False.

Returns

Weighted adjacency matrix.

Return type

numpy.ndarray

NEMtropy.matrix_generator.random_weighted_matrix_generator_uniform_custom_density(n, p=0.1, sup_ext=10, sym=False, seed=None, intweights=False)[source]

Generates a random weighted adjacency matrix with density determined by the parameter p and weights extracted from a uniform distribution.

Parameters
  • n (int) – Number of nodes.

  • p (float, optional) – Probability of observing a link between a couple of nodes, defaults to 0.1.

  • sup_ext (int or float, optional) – Maximum strength, defaults to 10.

  • sym (bool, optional) – If True the adjacency matrix is simmetric, defaults to False.

  • seed (int, optional) – Seed of the random process, defaults to None.

  • intweights (bool, optional) – If True weights are integers, defaults to False.

Returns

Weighted adjacency matrix.

Return type

numpy.ndarray

Model specific functions

NEMtropy.models_functions.expected_decm(x)
Expected parameters after the DECM.

It returns a concatenated array of out-degrees, in-degrees, out-strengths, in-strengths.

Parameters

x (numpy.ndarray) – DECM solution.

Returns

Expected parameters sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_decm_exp(theta)

Expected parameters after the DBCM. It returns a concatenated array of out-degrees and in-degrees. It is based on DBCM exponential version.

Parameters

theta – DBCM solution.

Returns

DBCM expected parameters sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_degree_cm(sol)

Computes the expected degrees of UBCM given the solution x.

Parameters

sol (numpy.ndarray) – UBCM solutions.

Returns

Expected degrees sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_ecm(sol)

Computes expected degrees and strengths sequence given solution x and y of UECM.

Parameters

sol (numpy.ndarray) – UECM solutions.

Returns

expected degrees and strengths sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_in_degree_dcm(sol)

Expected in-degree after the DBCM.

Parameters

sol (numpy.ndarray) – DBCM solution.

Returns

Expected in-degree sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_in_degree_dcm_exp(theta)

Expected in-degrees after the DBCM. It is based on DBCM exponential version.

Parameters

sol (numpy.ndarray) – DBCM solution.

Returns

In-degrees DBCM expectation.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_in_stregth_crema_directed(sol, adj)

Expected in-strength after CReMa.

Parameters
  • sol (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Tuple containing the original topology edges list and link related weigths.

Returns

Expected in-strength sequence

Return type

numpy.ndarray

NEMtropy.models_functions.expected_in_stregth_crema_directed_sparse(sol, adj)

Expected in-strength after CReMa. Sparse inisialization version.

Parameters
  • sol (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Tuple containing the original topology edges list and link related weigths.

Returns

Expected in-strength sequence

Return type

numpy.ndarray

NEMtropy.models_functions.expected_out_degree_dcm(sol)

Expected out-degree after the DBCM.

Parameters

sol (numpy.ndarray) – DBCM solution.

Returns

Expected out-degree sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_out_degree_dcm_exp(sol)

Expected out-degrees after the DBCM. It is based on DBCM exponential version.

Parameters

sol (numpy.ndarray) – DBCM solution.

Returns

Out-degrees DBCM expectation.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_out_strength_crema_directed(sol, adj)

Expected out-strength after CReMa.

Parameters
  • sol (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Tuple containing the original topology edges list and link related weigths.

Returns

Expected out-strength sequence

Return type

numpy.ndarray

NEMtropy.models_functions.expected_out_strength_crema_directed_sparse(sol, adj)

Expected out-strength after CReMa. Sparse initialisation version.

Parameters
  • sol (numpy.ndarray) – CReMa solution

  • adj ((numpy.ndarray, numpy.ndarray)) – Tuple containing the binary problem solution.

Returns

Expected out-strength sequence

Return type

numpy.ndarray

NEMtropy.models_functions.expected_strength_crema_undirected(sol, adj)

Computes the expected strengths of CReMa given its solution beta.

Parameters
  • sol (numpy.ndarray) – CReMa solutions.

  • adj (numpy.ndarray) – adjacency/pmatrix.

Returns

expected strengths sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.expected_strength_crema_undirected_sparse(sol, adj)

Computes the expected strengths of CReMa given its solution beta and the solutions of UBCM.

Parameters
  • sol (numpy.ndarray) – CReMa solutions.

  • adj (numpy.ndarray) – UBCM solutions.

Returns

expected strengths sequence.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_cm(x, args)

Returns the next UBCM iterative step for the fixed-point method.

Parameters
  • x (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray)) – Degrees and classes cardinality sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_cm_exp(x, args)

Returns the next UBCM iterative step for the fixed-point method. It is based on the exponential version of the UBCM.

Parameters
  • x (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray)) – Degrees and classes cardinality sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_crema_directed(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, adjacency binary/probability matrix, and non zero out and in indices

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_crema_directed_sparse(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The DBCM pmatrix is computed inside the function. Alternative version not in use.

Parameters
  • beta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Next iterative step.

Return type

numpy.ndarray]

NEMtropy.models_functions.iterative_crema_directed_sparse_2(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The DBCM pmatrix is computed inside the function. Alternative version not in use.

Parameters
  • beta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Next iterative step.

Return type

numpy.ndarray]

NEMtropy.models_functions.iterative_crema_sparse_2(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The UBCM pmatrix is computed inside the function. Alternative version not in use.

Parameters
  • beta (numpy.ndarray) – Previous iterative step..

  • args ((numpy.ndarray, numpy.ndarray)) – Strengths sequence and adjacency matrix.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_crema_undirected(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The UBCM pmatrix is pre-compute and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray)) – Strengths sequence and adjacency binary/probability matrix.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_crema_undirected_sparse(beta, args)

Returns the next CReMa iterative step for the fixed-point method. The UBCM pmatrix is computed inside the function.

Parameters
  • beta (numpy.ndarray) – Previous iterative step..

  • args ((numpy.ndarray, numpy.ndarray)) – Strengths sequence and adjacency matrix.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_dcm(x, args)

Returns the next DBCM iterative step for the fixed-point method.

Parameters
  • x (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, non zero out and in indices, and classes cardinalities sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_dcm_exp(theta, args)
Returns the next DBCM iterative step for the fixed-point 1 2.

It is based on the exponential version of the DBCM. This version only runs on non-zero indices.

Parameters
  • theta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Next iterative step.

Return type

numpy.ndarray

1

Squartini, Tiziano, and Diego Garlaschelli. “Analytical maximum-likelihood method to detect patterns in real networks.” New Journal of Physics 13.8 (2011): 083001. https://arxiv.org/abs/1103.0701

2

Squartini, Tiziano, Rossana Mastrandrea, and Diego Garlaschelli. “Unbiased sampling of network ensembles.” New Journal of Physics 17.2 (2015): 023052. https://arxiv.org/abs/1406.1197

NEMtropy.models_functions.iterative_dcm_exp_2(theta, args)
Returns the next DBCM iterative step for the fixed-point.

It is based on the exponential version of the DBCM. This version only runs on non-zero indices.

Parameters
  • theta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_decm(x, args)

Returns the next iterative step for the DECM Model.

Param

Previous iterative step.

Type

numpy.ndarray

Parameters

args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in degrees sequences, and out and in strengths sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_decm_exp(theta, args)

Returns the next iterative step for the DECM. It is based on the exponential version of the DBCM.

Parameters
  • theta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in degrees sequences, and out and in strengths sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_decm_exp_2(theta, args)

Returns the next iterative step for the DECM. It is based on the exponential version of the DBCM.

Parameters
  • theta (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Out and in degrees sequences, and out and in strengths sequences..

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_ecm(sol, args)

Returns the next UECM iterative step for the fixed-point method.

Parameters
  • sol (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray)) – Degrees and strengths sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.iterative_ecm_exp(sol, args)

Returns the next DECM iterative step for the fixed-point method. It is based on the exponential version of the UECM.

Parameters
  • sol (numpy.ndarray) – Previous iterative step.

  • args ((numpy.ndarray, numpy.ndarray)) – Degrees and strengths sequences.

Returns

Next iterative step.

Return type

numpy.ndarray

NEMtropy.models_functions.linsearch_fun_CM(X, args)

Linsearch function for UBCM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_CM_exp(X, args)

Linsearch function for UBCM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on UBCM exponential version.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_CM_exp_fixed(X)

Linsearch function for UBCM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on UBCM exponential version.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_CM_fixed(X)

Linsearch function for UBCM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DCM(xx, args)

Linsearch function for DBCM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DCM_exp(X, args)

Linsearch function for DBCM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on DBCM exponential version.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DCM_exp_fixed(X)

Linsearch function for DBCM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on DBCM exponential version.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DCM_fixed(xx)

Linsearch function for DBCM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

xx ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DECM(xx, args)

Linsearch function for DECM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DECM_exp(X, args)

Linsearch function for DECM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on DECM exponential version.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DECM_exp_fixed(X)

Linsearch function for DECM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on DECM exponential version.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_DECM_fixed(xx)

Linsearch function for DECM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_ECM(X, args)

Linsearch function for UECM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_ECM_exp(X, args)

Linsearch function for UECM newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on UBCM exponential version.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_ECM_exp_fixed(X)

Linsearch function for UECM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm. This function works on UBCM exponential version.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_ECM_fixed(X)

Linsearch function for UECM fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_crema_directed(xx, args)

Linsearch function for CReMa newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • xx ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f.

  • args ((function, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_crema_directed_fixed(xx)

Linsearch function for CReMa fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

xx ((numpy.ndarray, numpy.ndarray, numpy.ndarray, float, float, float)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_crema_undirected(X, args)

Linsearch function for CReMa newton and quasinewton methods. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters
  • X ((numpy.ndarray, numpy.ndarray, float, float, func)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, function f

  • args ((func, tuple)) – Tuple, step function and arguments.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.linsearch_fun_crema_undirected_fixed(X)

Linsearch function for CReMa fixed-point method. The function returns the step’s size, alpha. Alpha determines how much to move on the descending direction found by the algorithm.

Parameters

X ((numpy.ndarray, numpy.ndarray, float, float, int)) – Tuple of arguments to find alpha: solution, solution step, tuning parameter beta, initial alpha, step.

Returns

Working alpha.

Return type

float

NEMtropy.models_functions.loglikelihood_cm(x, args)

Returns UBCM loglikelihood function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Degrees and classes cardinality sequences.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_cm_exp(x, args)

Returns UBCM loglikelihood function evaluated in x. It is based on the exponential version of the UBCM.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Degrees and classes cardinality sequences.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_crema_directed(beta, args)

Returns CReMa loglikelihood function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_crema_directed_sparse(beta, args)

Returns CReMa loglikelihood function evaluated in beta. The DBCM pmatrix is computed inside the function. Sparse initialisation version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_crema_undirected(beta, args)

Returns CReMa loglikelihood function evaluated in beta. The UBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_crema_undirected_sparse(beta, args)

Computes CReMa loglikelihood function evaluated in beta. The UBCM pmatrix is computed inside the function. Sparse initialisation version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Strengths sequence and adjacency matrix.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_dcm(x, args)

Returns DBCM loglikelihood function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, non zero out and in indices, and classes cardinality sequence.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_dcm_exp(theta, args)

Returns DBCM * loglikelihood function evaluated in theta. It is based on the exponential version of the DBCM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and non zero out and in indices, and classes cardinalities sequence.

Returns

Loglikelihood value.

Return type

float

*

Squartini, Tiziano, and Diego Garlaschelli. “Analytical maximum-likelihood method to detect patterns in real networks.” New Journal of Physics 13.8 (2011): 083001. https://arxiv.org/abs/1103.0701

Squartini, Tiziano, Rossana Mastrandrea, and Diego Garlaschelli. “Unbiased sampling of network ensembles.” New Journal of Physics 17.2 (2015): 023052. https://arxiv.org/abs/1406.1197

NEMtropy.models_functions.loglikelihood_decm(x, args)

Returns DECM loglikelihood function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_decm_exp(x, args)

Returns DECM loglikelihood function evaluated in theta. It is based on the exponential version of the DECM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and out and in strengths sequences

Returns

Loglikelihood value.

Return type

float

Parisi, Federica, Tiziano Squartini, and Diego Garlaschelli. “A faster horse on a safer trail: generalized inference for the efficient reconstruction of weighted networks.” New Journal of Physics 22.5 (2020): 053053. https://arxiv.org/abs/1811.09829

NEMtropy.models_functions.loglikelihood_ecm(sol, args)

Returns UECM loglikelihood function evaluated in sol.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Degrees and strengths sequences.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_ecm_exp(sol, args)

Returns UECM loglikelihood function evaluated in sol. It is based on the exponential version of the UECM.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Degrees and strengths sequences.

Returns

Loglikelihood value.

Return type

float

NEMtropy.models_functions.loglikelihood_hessian_cm(x, args)

Returns UBCM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and classes cardinality sequences.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_cm_exp(x, args)

Returns UBCM loglikelihood hessian function evaluated in beta. It is based on the exponential version of the UBCM.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and classes cardinality sequences.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_crema_directed(beta, args)

Returns CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_crema_undirected(beta, args)

Returns CReMa loglikelihood hessian function evaluated in beta. The UBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_dcm(x, args)

Returns DBCM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Tuple containing out and in degrees sequences, non zero out and in indices, and classes cardinality sequence.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_dcm_exp(theta, args)

Returns DBCM § loglikelihood hessian function evaluated in theta. It is based on the exponential version of the DBCM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and non zero out and in indices, and the sequence of classes cardinalities.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

§

Squartini, Tiziano, and Diego Garlaschelli. “Analytical maximum-likelihood method to detect patterns in real networks.” New Journal of Physics 13.8 (2011): 083001. https://arxiv.org/abs/1103.0701

Squartini, Tiziano, Rossana Mastrandrea, and Diego Garlaschelli. “Unbiased sampling of network ensembles.” New Journal of Physics 17.2 (2015): 023052. https://arxiv.org/abs/1406.1197

NEMtropy.models_functions.loglikelihood_hessian_decm(x, args)

Returns DECM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_decm_exp(theta, args)

Returns DECM # loglikelihood hessian function evaluated in theta. It is based on the exponential version of the DECM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and out and in strengths sequences..

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

#

Parisi, Federica, Tiziano Squartini, and Diego Garlaschelli. “A faster horse on a safer trail: generalized inference for the efficient reconstruction of weighted networks.” New Journal of Physics 22.5 (2020): 053053. https://arxiv.org/abs/1811.09829

NEMtropy.models_functions.loglikelihood_hessian_diag_cm(x, args)

Returns the diagonal of the UBCM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and classes cardinality sequences.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_cm_exp(x, args)

Returns the diagonal of the UBCM loglikelihood hessian function evaluated in x. It is based on the exponential version of the UBCM.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and classes cardinality sequences.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_directed(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_directed_sparse(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is computed inside the function. Sparse initialization version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_directed_sparse_2(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed. Sparse initialization version. Alternative version not in use.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_sparse_2(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The UBCM pmatrix is computed inside the function. Sparse initialization version. Alternative version not in use.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_undirected(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_crema_undirected_sparse(beta, args)

Returns the diagonal of CReMa loglikelihood hessian function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed. Sparse initialization version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_dcm(x, args)

Returns the diagonal of DBCM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Out and in degrees sequences, non zero out and in indices, and classes cardinality sequence.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_dcm_exp(theta, args)

Returns the diagonal of the DBCM loglikelihood hessian function evaluated in theta. It is based on DBCM exponential version.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and non zero out and in indices, and the sequence of classes cardinalities.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

Squartini, Tiziano, and Diego Garlaschelli. “Analytical maximum-likelihood method to detect patterns in real networks.” New Journal of Physics 13.8 (2011): 083001. https://arxiv.org/abs/1103.0701

Squartini, Tiziano, Rossana Mastrandrea, and Diego Garlaschelli. “Unbiased sampling of network ensembles.” New Journal of Physics 17.2 (2015): 023052. https://arxiv.org/abs/1406.1197

NEMtropy.models_functions.loglikelihood_hessian_diag_decm(x, args)

Returns the diagonal of DECM loglikelihood hessian function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_decm_exp(theta, args)

Returns the diagonal of the DECM loglikelihood hessian function evaluated in theta. It is based on the DECM exponential version.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood hessian diagonal.

Return type

numpy.ndarray

Parisi, Federica, Tiziano Squartini, and Diego Garlaschelli. “A faster horse on a safer trail: generalized inference for the efficient reconstruction of weighted networks.” New Journal of Physics 22.5 (2020): 053053. https://arxiv.org/abs/1811.09829

NEMtropy.models_functions.loglikelihood_hessian_diag_ecm(sol, args)

Returns the diagonal of UECM loglikelihood hessian function evaluated in sol.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and strengths sequences.

Returns

Hessian matrix diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_diag_ecm_exp(sol, args)

Returns the diagonal of UECM loglikelihood hessian function evaluated in sol. It is based on the exponential version of the UECM.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian function. Degrees and strengths sequences.

Returns

Hessian matrix diagonal.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_ecm(sol, args)

Returns DBCM loglikelihood hessian function evaluated in sol.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian. Degrees and strengths sequences.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_hessian_ecm_exp(sol, args)

Returns DBCM loglikelihood hessian function evaluated in sol. It is based on the exponential version of the UECM.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood hessian. Degrees and strengths sequences.

Returns

Loglikelihood hessian matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_cm(x, args)

Returns UBCM loglikelihood gradient function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Degrees and classes cardinality sequences.

Returns

Loglikelihood gradient.

Return type

float

NEMtropy.models_functions.loglikelihood_prime_cm_exp(x, args)

Returns UBCM loglikelihood gradient function evaluated in beta. It is based on the exponential version of the UBCM.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Degrees and classes cardinality sequences.

Returns

Loglikelihood gradient.

Return type

float

NEMtropy.models_functions.loglikelihood_prime_crema_directed(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood gradient value.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_crema_directed_sparse(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The DBCM pmatrix is pre-computed and explicitly passed. Sparse initialization version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_crema_directed_sparse_2(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The DBCM pmatrix is computed inside the function. Sparse initialization version. Alternative version not used.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Out and in strengths sequences, adjacency matrix, and non zero out and in indices.

Returns

Loglikelihood gradient value.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_crema_sparse_2(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The UBCM pmatrix is computed inside the function. Sparse initialization version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood gradient value.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_crema_undirected(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The UBCM pmatrix is pre-computed and explicitly passed.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood gradient value.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_crema_undirected_sparse(beta, args)

Returns CReMa loglikelihood gradient function evaluated in beta. The UBCM pmatrix is pre-computed and explicitly passed. Sparse initialization version.

Parameters
  • beta (numpy.ndarray) – Evaluating point beta.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient function. Strengths sequence and adjacency binary/probability matrix.

Returns

Loglikelihood gradient value.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_dcm(x, args)

Returns DBCM loglikelihood gradient function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient. Out and in degrees sequences, non zero out and in indices, and classes cardinality sequence.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_dcm_exp(theta, args)

Returns DBCM ** loglikelihood gradient function evaluated in theta. It is based on the exponential version of the DBCM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient. Out and in degrees sequences, and non zero out and in indices, and the sequence of classes cardinalities.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

Squartini, Tiziano, and Diego Garlaschelli. “Analytical maximum-likelihood method to detect patterns in real networks.” New Journal of Physics 13.8 (2011): 083001. https://arxiv.org/abs/1103.0701

**

Squartini, Tiziano, Rossana Mastrandrea, and Diego Garlaschelli. “Unbiased sampling of network ensembles.” New Journal of Physics 17.2 (2015): 023052. https://arxiv.org/abs/1406.1197

NEMtropy.models_functions.loglikelihood_prime_decm(x, args)

Returns DECM loglikelihood gradient function evaluated in x.

Parameters
  • x (numpy.ndarray) – Evaluating point x.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood gradient.

Return type

numpy.array

NEMtropy.models_functions.loglikelihood_prime_decm_exp(theta, args)[source]

Returns DECM †† loglikelihood gradient function evaluated in theta. It is based on the exponential version of the DECM.

Parameters
  • theta (numpy.ndarray) – Evaluating point theta.

  • args ((numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood function. Out and in degrees sequences, and out and in strengths sequences.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

††

Parisi, Federica, Tiziano Squartini, and Diego Garlaschelli. “A faster horse on a safer trail: generalized inference for the efficient reconstruction of weighted networks.” New Journal of Physics 22.5 (2020): 053053. https://arxiv.org/abs/1811.09829

NEMtropy.models_functions.loglikelihood_prime_ecm(sol, args)

Returns DECM loglikelihood gradient function evaluated in sol.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient. Degrees and strengths sequences.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

NEMtropy.models_functions.loglikelihood_prime_ecm_exp(sol, args)

Returns DECM loglikelihood gradient function evaluated in sol. It is based on the exponential version of the UECM.

Parameters
  • sol (numpy.ndarray) – Evaluating point sol.

  • args ((numpy.ndarray, numpy.ndarray)) – Arguments to define the loglikelihood gradient. Degrees and strengths sequences.

Returns

Loglikelihood gradient.

Return type

numpy.ndarray

NEMtropy.models_functions.pmatrix_cm(x, args)[source]

Computes and returns the probability matrix induced by UBCM.

Parameters
  • x (numpy.ndarray) – Solutions of UBCM.

  • args ((int, )) – Number of nodes.

Returns

UBCM probability matrix.

Return type

numpy.ndarray

NEMtropy.models_functions.pmatrix_dcm(x, args)

Computes and returns the probability matrix induced by DBCM.

Parameters
  • x (numpy.ndarray) – DBCM solution.

  • args ((int, numpy.ndarray, numpy.ndarray)) – Problem dimension, out and in indices of non-zero nodes.

Returns

DBCM probability matrix

Return type

numpy.ndarray

Solver functions

NEMtropy.solver_functions.matrix_regulariser_function(b, eps)[source]

Trasforms input matrix in a positive defined matrix by adding positive quantites to the main diagonal.

Parameters
  • b (numpy.ndarray) – Matrix.

  • eps (float) – Positive quantity to add.

Returns

Regularised matrix.

Return type

numpy.ndarray

NEMtropy.solver_functions.matrix_regulariser_function_eigen_based(b, eps)[source]

Trasform input matrix in a positive defined matrix by regularising eigenvalues.

Parameters
  • b (numpy.ndarray) – Matrix.

  • eps (float) – Positive quantity to add.

Returns

Regularised matrix.

Return type

numpy.ndarray

NEMtropy.solver_functions.solver(x0, fun, step_fun, linsearch_fun, hessian_regulariser, fun_jac=None, tol=1e-06, eps=1e-10, max_steps=100, method='newton', verbose=False, regularise=True, regularise_eps=0.001, full_return=False, linsearch=True)[source]

Find roots of eq. fun = 0, using newton, quasinewton or fixed-point algorithm.

Parameters
  • x0 (numpy.ndarray) – Initial point

  • fun (function) – Function handle of the function to find the roots of.

  • step_fun (function) – Function to compute the algorithm step

  • linsearch_fun (function) – Function to compute the linsearch

  • hessian_regulariser (function) – Function to regularise fun hessian

  • fun_jac (function, optional) – Function to compute the hessian of fun, defaults to None

  • tol (float, optional) – The solver stops when |fun|<tol, defaults to 1e-6

  • eps (float, optional) – The solver stops when the difference between two consecutive steps is less than eps, defaults to 1e-10

  • max_steps (int, optional) – Maximum number of steps the solver takes, defaults to 100

  • method (str, optional) – Method the solver uses to solve the problem. Choices are “newton”, “quasinewton”, “fixed-point”. Defaults to “newton”

  • verbose (bool, optional) – If True the solver prints out information at each step, defaults to False

  • regularise (bool, optional) – If True the solver will regularise the hessian matrix, defaults to True

  • regularise_eps (float, optional) – Positive value to pass to the regulariser function, defaults to 1e-3

  • full_return (bool, optional) – If True the function returns information on the convergence, defaults to False

  • linsearch (bool, optional) – If True a linsearch algorithm is implemented, defaults to True

Returns

Solution to the optimization problem

Return type

numpy.ndarray

NEMtropy.solver_functions.sufficient_decrease_condition(f_old, f_new, alpha, grad_f, p, c1=0.0001, c2=0.9)

Returns True if upper wolfe condition is respected.

Parameters
  • f_old (float) – Function value at previous iteration.

  • f_new (float) – Function value at current iteration.

  • alpha (float) – Alpha parameter of linsearch.

  • grad_f (numpy.ndarray) – Function gradient.

  • p (numpy.ndarray) – Current iteration increment.

  • c1 (float, optional) – Tuning parameter, defaults to 1e-04.

  • c2 (float, optional) – Tuning parameter, defaults to 0.9.

Returns

Condition validity.

Return type

bool

Network functions

NEMtropy.network_functions.build_graph_fast(edgelist, is_directed)

Generates adjacency matrix given edgelist, numpy array format is used.

Parameters
  • edgelist (numpy.ndarray) – Edgelist.

  • is_directed (bool) – True if edge direction is informative.

Returns

Adjacency matrix.

Return type

numpy.ndarray

NEMtropy.network_functions.build_graph_fast_weighted(edgelist, is_directed)

Generates weighted adjacency matrix given edgelist, numpy array format is used.

Parameters
  • edgelist (numpy.ndarray) – Edgelist.

  • is_directed (bool) – True if edge direction is informative.

Returns

Adjacency matrix.

Return type

numpy.ndarray

NEMtropy.network_functions.build_graph_from_edgelist(edgelist, is_directed, is_sparse=False, is_weighted=False)[source]

Generates adjacency matrix given edgelist.

Parameters
  • edgelist (numpy.ndarray) – Edgelist.

  • is_directed (bool) – True if edge direction is informative.

  • is_sparse (bool) – If true the output is a sparse matrix.

  • is_weighted (bool) – True if the edgelist is weighted

Returns

Adjacency matrix.

Return type

numpy.ndarray or scipy.sparse.csr_matrix

NEMtropy.network_functions.build_graph_sparse(edgelist, is_directed)[source]

Generates adjacency matrix given edgelist, scipy sparse format is used.

Parameters
  • edgelist (numpy.ndarray) – Edgelist.

  • is_sparse (bool) – True if edge direction is informative.

Returns

Adjacency matrix.

Return type

scipy.sparse.csr_matrix

NEMtropy.network_functions.build_graph_sparse_weighted(edgelist, is_directed)[source]

Generates weighted adjacency matrix given edgelist, scipy sparse format is used.

Parameters
  • edgelist (numpy.ndarray) – Edgelist.

  • is_sparse (bool) – True if edge direction is informative.

Returns

Adjacency matrix.

Return type

scipy.sparse.csr_matrix

NEMtropy.network_functions.degree(a)[source]

Returns matrix a degrees sequence.

Parameters

a (numpy.ndarray, scipy.sparse) – Adjacency matrix.

Returns

Degree sequence.

Return type

numpy.ndarray

NEMtropy.network_functions.edgelist_from_edgelist_directed(edgelist)[source]

Creates a new edgelists replacing nodes labels with indexes. Returns also two dictionaries that keep track of the nodes index-label relation. Works also on weighted graphs.

Parameters

edgelist (list) – List of edges.

Returns

Re-indexed list of edges, out-degrees, in-degrees, index to label dictionary

Return type

(dict, numpy.ndarray, numpy.ndarray, dict)

NEMtropy.network_functions.edgelist_from_edgelist_undirected(edgelist)[source]

Creates a new edgelist with the indexes of the nodes instead of the names. Returns also a dictionary that keep track of the nodes and, depending on the type of graph, degree and strengths sequences.

Parameters

edgelist (numpy.ndarray or list) – edgelist.

Returns

edgelist, degrees sequence, strengths sequence and new labels to old labels dictionary.

Return type

(numpy.ndarray, numpy.ndarray, numpy.ndarray, dict)

NEMtropy.network_functions.in_degree(a)[source]

Returns matrix a in degrees sequence.

Parameters

a (numpy.ndarray, scipy.sparse.csr.csr_matrix, scipy.sparse.coo.coo_matrix) – Adjacency matrix.

Returns

In degree sequence.

Return type

numpy.ndarray

NEMtropy.network_functions.in_strength(a)[source]

Returns matrix a in strengths sequence.

Parameters

a (numpy.ndarray, scipy.sparse.csr.csr_matrix, scipy.sparse.coo.coo_matrix) – Adjacency matrix.

Returns

In strengths sequence.

Return type

numpy.ndarray

NEMtropy.network_functions.out_degree(a)[source]

Returns matrix a out degrees sequence.

Parameters

a (numpy.ndarray, scipy.sparse.csr.csr_matrix, scipy.sparse.coo.coo_matrix) – Adjacency matrix

Returns

Out degree sequence.

Return type

numpy.ndarray

NEMtropy.network_functions.out_strength(a)[source]

Returns matrix a out strengths sequence.

Parameters

a (numpy.ndarray, scipy.sparse.csr.csr_matrix, scipy.sparse.coo.coo_matrix) – Adjacency matrix.

Returns

Out strengths sequence.

Return type

numpy.ndarray

NEMtropy.network_functions.strength(a)[source]

Returns matrix a strengths sequence.

Parameters

a (numpy.ndarray, scipy.sparse) – Adjacency matrix.

Returns

Strengths sequence.

Return type

numpy.ndarray

Poibin class