Skip to content

mnns.dynamics

Functions to evolve nanowire networks over time.

Functions:

Name Description
get_evolution_current

Calculates the current draw from each drain node as a function of time,

get_evolution_node_voltages

To be used in conjunction with solve_evolution. Takes the output from

set_state_variables

Deprecated. Please use NanowireNetwork.set_state_var() instead.

solve_evolution

Deprecated. Please use NanowireNetwork.evolve() instead.

get_evolution_current

get_evolution_current(
    NWN: NanowireNetwork,
    sol: OdeResult,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    state_vars: Optional[list[str]] = None,
    scaled: bool = False,
    solver: str = "spsolve",
    **kwargs
) -> npt.NDArray

Calculates the current draw from each drain node as a function of time, given the output state variables from NWN.evolve().

The source, drain, and voltage are expected to be the same as in the evolution of the state variables.

This assumes the resistance function only depends on the first state variable set.

Parameters:

Name Type Description Default
NWN NanowireNetwork

Nanowire network.

required
sol OdeResult

Output state variables from NWN.evolve().

required
source_node NWNNode, or list of NWNNodes

Voltage source nodes.

required
drain_node NWNNode, or list of NWNNodes

Grounded output nodes.

required
voltage_func Callable

Function which inputs the time as a scalar and returns the voltage of all the source nodes as a scalar. Should be dimensionless.

required
state_vars list of str

List of state variables to use. Should not need to be provided. Defaults to NWN.state_vars.

None
scaled bool

Scale the output by the characteristic values. Default: False.

False
solver str

Name of sparse matrix solving algorithm to use. Default: "spsolve".

'spsolve'
**kwargs

Keyword arguments passed to the solver.

{}

Returns:

Name Type Description
current_array ndarray

Array containing the current flow through each drain node. Each column corresponds to a drain node in the order passed.

Source code in mnns/dynamics.py
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
def get_evolution_current(
    NWN: NanowireNetwork,
    sol: OdeResult,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    state_vars: Optional[list[str]] = None,
    scaled: bool = False,
    solver: str = "spsolve",
    **kwargs,
) -> npt.NDArray:
    """
    Calculates the current draw from each drain node as a function of time,
    given the output state variables from `NWN.evolve()`.

    The source, drain, and voltage are expected to be the same as in the
    evolution of the state variables.

    This assumes the resistance function only depends on the first state variable
    set.

    Parameters
    ----------
    NWN : NanowireNetwork
        Nanowire network.

    sol : OdeResult
        Output state variables from `NWN.evolve()`.

    source_node : NWNNode, or list of NWNNodes
        Voltage source nodes.

    drain_node : NWNNode, or list of NWNNodes
        Grounded output nodes.

    voltage_func : Callable
        Function which inputs the time as a scalar and returns the voltage of
        all the source nodes as a scalar. Should be dimensionless.

    state_vars : list of str, optional
        List of state variables to use. Should not need to be provided.
        Defaults to `NWN.state_vars`.

    scaled : bool, optional
        Scale the output by the characteristic values. Default: False.

    solver : str, optional
        Name of sparse matrix solving algorithm to use. Default: "spsolve".

    **kwargs
        Keyword arguments passed to the solver.

    Returns
    -------
    current_array: ndarray
        Array containing the current flow through each drain node. Each column
        corresponds to a drain node in the order passed.

    """
    # Get lists of source and drain nodes
    if isinstance(source_node, tuple):
        source_node = [source_node]
    if isinstance(drain_node, tuple):
        drain_node = [drain_node]

    if state_vars is None:
        state_vars = NWN.state_vars

    # Preallocate output
    current_array = np.zeros((len(sol.t), len(drain_node)))

    # Loop through each time step
    for i, state in enumerate(sol.y.T):

        # Set state variables
        split = np.split(state, len(state_vars))
        for j, var in enumerate(state_vars):
            NWN.set_state_var(var, split[j])

        # Update resistance
        NWN.update_resistance(split[0])

        # Find current draw through drain nodes
        V = voltage_func(sol.t[i])
        current_array[i] = solve_drain_current(
            NWN, source_node, drain_node, V, scaled, solver, **kwargs
        )

    return current_array.squeeze()

get_evolution_node_voltages

get_evolution_node_voltages(
    NWN: Graph,
    sol: OdeResult,
    edge_list: list[NWNEdge],
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    solver: str = "spsolve",
    **kwargs
) -> npt.NDArray

To be used in conjunction with solve_evolution. Takes the output from solve_evolution and finds the voltage of all nodes in the network at each time step.

The appropriate parameters passed should be the same as solve_evolution.

Parameters:

Name Type Description Default
NWN Graph

Nanowire network.

required
sol OdeResult

Output from solve_evolution.

required
edge_list list of tuples

Output from solve_evolution.

required
source_node tuple, or list of tuples

Voltage source nodes.

required
drain_node tuple, or list of tuples

Grounded output nodes.

required
voltage_func Callable

The applied voltage with the calling signature func(t). The voltage should have units of v0.

required
solver str

Name of sparse matrix solving algorithm to use. Default: "spsolve".

'spsolve'
**kwargs

Keyword arguments passed to the solver.

{}

Returns:

Name Type Description
node_voltage_array ndarray

An n x m array containing the voltage of each node in the network, corresponding to the n time steps and m nodes in the network.

Source code in mnns/dynamics.py
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
def get_evolution_node_voltages(
    NWN: nx.Graph,
    sol: OdeResult,
    edge_list: list[NWNEdge],
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    solver: str = "spsolve",
    **kwargs,
) -> npt.NDArray:
    """
    To be used in conjunction with `solve_evolution`. Takes the output from
    `solve_evolution` and finds the voltage of all nodes in the network at each
    time step.

    The appropriate parameters passed should be the same as `solve_evolution`.

    Parameters
    ----------
    NWN : Graph
        Nanowire network.

    sol : OdeResult
        Output from `solve_evolution`.

    edge_list : list of tuples
        Output from `solve_evolution`.

    source_node : tuple, or list of tuples
        Voltage source nodes.

    drain_node : tuple, or list of tuples
        Grounded output nodes.

    voltage_func : Callable
        The applied voltage with the calling signature `func(t)`. The voltage
        should have units of `v0`.

    solver : str, optional
        Name of sparse matrix solving algorithm to use. Default: "spsolve".

    **kwargs
        Keyword arguments passed to the solver.

    Returns
    -------
    node_voltage_array: ndarray
        An `n` x `m` array containing the voltage of each node in the network,
        corresponding to the `n` time steps and `m` nodes in the network.

    """
    # Get lists of source and drain nodes
    if isinstance(source_node, tuple):
        source_node = [source_node]
    if isinstance(drain_node, tuple):
        drain_node = [drain_node]

    # Preallocate output
    node_voltage_array = np.zeros((len(sol.t), len(NWN.nodes)))

    # Loop through each time step
    for i in range(len(sol.t)):
        # Set state variables and get drain currents
        input_V = voltage_func(sol.t[i])
        set_state_variables(NWN, sol.y.T[i], edge_list)
        *node_voltage_array[i], I = solve_network(
            NWN, source_node, drain_node, input_V, "voltage", solver, **kwargs
        )

    return node_voltage_array.squeeze()

set_state_variables

set_state_variables(NWN: NanowireNetwork, *args)

Deprecated. Please use NanowireNetwork.set_state_var() instead.

Sets the given nanowire network's state variable. Can be called in the following ways:

set_state_variables(NWN, w)
    where `w` is a scalar value which is set for all edges.

set_state_variables(NWN, w, edge_list)
    where `w` is an ndarray and edge_list is a list. The `w` array
    contains the state variable for the corresponding edge in
    `edge_list`.

set_state_variables(NWN, w, tau)
    where `w` and `tau` are scalar values which is set for all edges.

set_state_variables(NWN, w, tau, edge_list)
    where `w` and `tau` are ndarrays and edge_list is a list. The `w`
    and `tau` arrays contains the state variables for the
    corresponding edge in `edge_list`.

set_state_variables(NWN, w, tau, epsilon)
    where `w`, `tau`, and `epsilon` are scalar values which is set for
    all edges.

set_state_variables(NWN, w, tau, epsilon, edge_list)
    where `w`,`tau`, and `epsilon` are ndarrays and edge_list is a list.
    The `w`,`tau`, and `epsilon` arrays contains the state variables
    for the corresponding edge in `edge_list`.

Parameters:

Name Type Description Default
NWN NanowireNetwork

Nanowire network.

required
*args

See above.

()
Source code in mnns/dynamics.py
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
def set_state_variables(NWN: NanowireNetwork, *args):
    """
    Deprecated. Please use NanowireNetwork.set_state_var() instead.

    Sets the given nanowire network's state variable. Can be called in the
    following ways:

        set_state_variables(NWN, w)
            where `w` is a scalar value which is set for all edges.

        set_state_variables(NWN, w, edge_list)
            where `w` is an ndarray and edge_list is a list. The `w` array
            contains the state variable for the corresponding edge in
            `edge_list`.

        set_state_variables(NWN, w, tau)
            where `w` and `tau` are scalar values which is set for all edges.

        set_state_variables(NWN, w, tau, edge_list)
            where `w` and `tau` are ndarrays and edge_list is a list. The `w`
            and `tau` arrays contains the state variables for the
            corresponding edge in `edge_list`.

        set_state_variables(NWN, w, tau, epsilon)
            where `w`, `tau`, and `epsilon` are scalar values which is set for
            all edges.

        set_state_variables(NWN, w, tau, epsilon, edge_list)
            where `w`,`tau`, and `epsilon` are ndarrays and edge_list is a list.
            The `w`,`tau`, and `epsilon` arrays contains the state variables
            for the corresponding edge in `edge_list`.

    Parameters
    ----------
    NWN: Graph
        Nanowire network.

    *args
        See above.

    """
    # Only scalar w passed
    if len(args) == 1:
        if isinstance(args[0], Number):
            w = args[0]
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {"w": w, "conductance": 1 / R}
                for edge in NWN.edges
                if NWN.edges[edge]["type"] == "junction"
            }
            nx.set_edge_attributes(NWN, attrs)

        else:
            raise ValueError("Invalid arguments.")

    elif len(args) == 2:
        # vector w and edge_list passed
        if isinstance(args[0], np.ndarray) and isinstance(args[1], Iterable):
            w, edge_list = args
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {"w": w[i], "conductance": 1 / R[i]}
                for i, edge in enumerate(edge_list)
            }
            nx.set_edge_attributes(NWN, attrs)

        # scalar w and scalar tau passed
        elif isinstance(args[0], Number) and isinstance(args[1], Number):
            w, tau = args
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {"w": w, "conductance": 1 / R, "tau": tau}
                for edge in NWN.edges
                if NWN.edges[edge]["type"] == "junction"
            }
            nx.set_edge_attributes(NWN, attrs)
            NWN.graph["tau"] = tau

        else:
            raise ValueError("Invalid arguments.")

    elif len(args) == 3:
        # vector w, vector tau, and edge_list passed
        if (
            isinstance(args[0], np.ndarray)
            and isinstance(args[1], np.ndarray)
            and isinstance(args[2], Iterable)
        ):
            w, tau, edge_list = args
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {"w": w[i], "conductance": 1 / R[i], "tau": tau[i]}
                for i, edge in enumerate(edge_list)
            }
            nx.set_edge_attributes(NWN, attrs)

        # scalar w, scalar tau, scalar epsilon passed
        if (
            isinstance(args[0], Number)
            and isinstance(args[1], Number)
            and isinstance(args[2], Number)
        ):
            w, tau, epsilon = args
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {
                    "w": w,
                    "conductance": 1 / R,
                    "tau": tau,
                    "epsilon": epsilon,
                }
                for edge in NWN.edges
                if NWN.edges[edge]["type"] == "junction"
            }
            nx.set_edge_attributes(NWN, attrs)
            NWN.graph["tau"] = tau
            NWN.graph["epsilon"] = epsilon

        else:
            raise ValueError("Invalid arguments.")

    elif len(args) == 4:
        # vector w, vector tau, vector epsilon, and edge_list passed
        if (
            isinstance(args[0], np.ndarray)
            and isinstance(args[1], np.ndarray)
            and isinstance(args[2], np.ndarray)
            and isinstance(args[3], Iterable)
        ):
            w, tau, epsilon, edge_list = args
            R = models.resist_func(NWN, w)

            attrs = {
                edge: {
                    "w": w[i],
                    "conductance": 1 / R[i],
                    "tau": tau[i],
                    "epsilon": epsilon[i],
                }
                for i, edge in enumerate(edge_list)
            }
            nx.set_edge_attributes(NWN, attrs)

        else:
            raise ValueError("Invalid arguments.")

    else:
        raise ValueError("Invalid number of arguments.")

solve_evolution

solve_evolution(
    NWN: NanowireNetwork,
    t_eval: NDArray,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable = None,
    tol: float = 1e-12,
    model: str = "default",
    solver: str = "spsolve",
    **kwargs
) -> tuple[OdeResult, list[NWNEdge]]

Deprecated. Please use NanowireNetwork.evolve() instead.

Solve for the state variables w of the junctions of the given nanowire network at various points in time with an applied voltage.

Parameters:

Name Type Description Default
NWN NanowireNetwork

Nanowire network.

required
t_eval ndarray

Time points to evaluate the nanowire network at. These should have units of t0.

required
source_node tuple, or list of tuples

Voltage source nodes.

required
drain_node tuple, or list of tuples

Grounded output nodes.

required
voltage_func Callable

The applied voltage with the calling signature func(t). The voltage should have units of v0.

required
window_func Callable

The window function used in the derivative of w. The calling signature is f(w) where w is the array of state variables. The default window function is f(w) = 1.

None
tol float

Tolerance of scipy.integrate.solve_ivp. Defaults to 1e-12.

1e-12
model (default, decay, chen)

Evolutionary model type. Default: "default".

"default"
solver str

Name of sparse matrix solving algorithm to use. Default: "spsolve".

'spsolve'
**kwargs

Keyword arguments passed to the solver.

{}

Returns:

Name Type Description
sol OdeResult

Output from scipy.integrate.solve_ivp. See the SciPy documentation for information on this output's formatting.

edge_list list of tuples

List of the edges corresponding with each w.

Source code in mnns/dynamics.py
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
def solve_evolution(
    NWN: NanowireNetwork,
    t_eval: npt.NDArray,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable = None,
    tol: float = 1e-12,
    model: str = "default",
    solver: str = "spsolve",
    **kwargs,
) -> tuple[OdeResult, list[NWNEdge]]:
    """
    Deprecated. Please use NanowireNetwork.evolve() instead.

    Solve for the state variables `w` of the junctions of the given nanowire
    network at various points in time with an applied voltage.

    Parameters
    ----------
    NWN: Graph
        Nanowire network.

    t_eval : ndarray
        Time points to evaluate the nanowire network at. These should have
        units of `t0`.

    source_node : tuple, or list of tuples
        Voltage source nodes.

    drain_node : tuple, or list of tuples
        Grounded output nodes.

    voltage_func : Callable
        The applied voltage with the calling signature `func(t)`. The voltage
        should have units of `v0`.

    window_func : Callable, optional
        The window function used in the derivative of `w`. The calling
        signature is `f(w)` where w is the array of state variables.
        The default window function is `f(w) = 1`.

    tol : float, optional
        Tolerance of `scipy.integrate.solve_ivp`. Defaults to 1e-12.

    model : {"default", "decay", "chen"}, optional
        Evolutionary model type. Default: "default".

    solver : str, optional
        Name of sparse matrix solving algorithm to use. Default: "spsolve".

    **kwargs
        Keyword arguments passed to the solver.

    Returns
    -------
    sol : OdeResult
        Output from `scipy.integrate.solve_ivp`. See the SciPy documentation
        for information on this output's formatting.

    edge_list : list of tuples
        List of the edges corresponding with each `w`.

    """
    if model not in ("default", "HP", "decay", "chen", "SLT"):
        raise ValueError(f'Unsupported model type: model="{model}"')

    # Default window function
    if window_func is None:
        window_func = lambda x: 1

    # Get list of junction edges and state variable w
    edge_list, w0 = map(
        list,
        zip(*[((u, v), w) for u, v, w in NWN.edges.data("w") if w is not None]),
    )

    # Get edge indices
    start_nodes, end_nodes = get_edge_indices(NWN, edge_list)

    # Get list of tau
    tau0 = [tau for _, _, tau in NWN.edges.data("tau") if tau is not None]

    # Get list of epsilon
    epsilon0 = [ep for _, _, ep in NWN.edges.data("epsilon") if ep is not None]

    # Define initial state and associated derivative
    if model == "default" or model == "HP":
        _deriv = models.HP_model
        y0 = w0
    elif model == "decay":
        _deriv = models.decay_HP_model
        y0 = w0
    elif model == "chen" or model == "SLT":
        _deriv = models.SLT_HP_model
        y0 = np.hstack((w0, tau0, epsilon0))
        if "sigma" not in NWN.graph.keys():
            raise AttributeError(
                "sigma, theta, and a must be set before using the Chen model."
            )

    # Solve the system of ODEs
    t_span = (t_eval[0], t_eval[-1])
    sol = solve_ivp(
        _deriv,
        t_span,
        y0,
        "DOP853",
        t_eval,
        atol=tol,
        rtol=tol,
        args=(
            NWN,
            source_node,
            drain_node,
            voltage_func,
            window_func,
            solver,
            kwargs,
        ),
    )

    # Update NWN to final state variables.
    if model == "default":
        set_state_variables(NWN, sol.y[:, -1], edge_list)
    elif model == "decay":
        set_state_variables(NWN, sol.y[:, -1], edge_list)
    elif model == "chen":
        w_list, tau_list, eps_list = np.split(sol.y, 3)
        set_state_variables(
            NWN, w_list[:, -1], tau_list[:, -1], eps_list[:, -1], edge_list
        )

    return sol, edge_list