Skip to content

mnns.models

Various dynamic models for nanowire networks.

Functions:

Name Description
HP_model

HP Model [1]. Provides the time derivative of the state variable x (the

SLT_HP_model

SLT HP Model [1]. Provides the time derivative of the state variable x

decay_HP_model

Decay HP Model [1]. Provides the time derivative of the state variable x

resist_func

Linear resistance function in nondimensionalized form.

HP_model

HP_model(
    t: float,
    x: NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: Optional[dict] = None,
) -> npt.NDArray

HP Model [1]. Provides the time derivative of the state variable x (the dimensionless version of w). Assumes voltage sources are used.

Parameters:

Name Type Description Default
t float

Current time to solve at.

required
x ndarray

Array containing the state variable x for each junction in the NWN.

required
NWN NanowireNetwork

Input nanowire network graph.

required
source_node NWNNode or list of NWNNode

Source node(s) of the input voltage.

required
drain_node NWNNode or list of NWNNode

Drain/grounded node(s).

required
voltage_func Callable

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

required
window_func Callable

Function which inputs the state variable x as an array and returns the window function value as an array.

required
solver str

SciPy sparse matrix equation solver.

'spsolve'
**kwargs Optional[dict]

Keyword arguments to pass to the SciPy sparse matrix equation solver.

None

Returns:

Name Type Description
dxdt ndarray

Array of the time derivative of the state variable x.

References

[1] D. B. Strukov, G. S. Snider, D. R. Stewart and R. S. Williams, Nature, 2008, 453, 80-83

Source code in mnns/models.py
 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
def HP_model(
    t: float,
    x: npt.NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: Optional[dict] = None,
) -> npt.NDArray:
    """
    HP Model [1]. Provides the time derivative of the state variable `x` (the
    dimensionless version of `w`). Assumes voltage sources are used.

    Parameters
    ----------
    t : float
        Current time to solve at.

    x : ndarray
        Array containing the state variable `x` for each junction in the NWN.

    NWN : NanowireNetwork
        Input nanowire network graph.

    source_node : NWNNode or list of NWNNode
        Source node(s) of the input voltage.

    drain_node : NWNNode or list of NWNNode
        Drain/grounded node(s).

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

    window_func : Callable
        Function which inputs the state variable `x` as an array and returns
        the window function value as an array.

    solver : str
        SciPy sparse matrix equation solver.

    **kwargs
        Keyword arguments to pass to the SciPy sparse matrix equation solver.

    Returns
    -------
    dxdt : ndarray
        Array of the time derivative of the state variable `x`.

    References
    ----------
    [1] D. B. Strukov, G. S. Snider, D. R. Stewart and R. S. Williams,
        *Nature*, 2008, **453**, 80-83

    """
    if kwargs is None:
        kwargs = dict()

    # Update all wire junction resistances
    R = NWN.update_resistance(x)

    # Find applied voltage at the current time
    applied_V = voltage_func(t)

    # Solve for voltage at each node
    *V, I = solve_network(
        NWN, source_node, drain_node, applied_V, "voltage", solver, **kwargs
    )
    V = np.asarray(V)

    # Get start and end indices
    start, end = NWN.wire_junction_indices()

    # Find voltage differences
    v0 = V[start]
    v1 = V[end]
    V_delta = np.abs(v0 - v1) * np.sign(applied_V)

    # Find dw/dt
    dxdt = V_delta / R * window_func(x)

    return dxdt

SLT_HP_model

SLT_HP_model(
    t: float,
    y: NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: Optional[dict] = None,
) -> npt.NDArray

SLT HP Model [1]. Provides the time derivative of the state variable x (the dimensionless version of w), tau, and epsilon. Assumes that these state variable are concatenated into a single 1D array input.

Assumes voltage sources are used.

Requires NWN.graph["tau"] to be set to the decay constant value.

Parameters:

Name Type Description Default
t float

Current time to solve at.

required
y ndarray

Array containing the state variables x, tau, and epsilon for each junction in the NWN concatenated into a single 1D array input.

required
NWN NanowireNetwork

Input nanowire network graph.

required
source_node NWNNode or list of NWNNode

Source node(s) of the input voltage.

required
drain_node NWNNode or list of NWNNode

Drain/grounded node(s).

required
voltage_func Callable

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

required
window_func Callable

Function which inputs the state variable x as an array and returns the window function value as an array.

required
solver str

SciPy sparse matrix equation solver.

'spsolve'
**kwargs Optional[dict]

Keyword arguments to pass to the SciPy sparse matrix equation solver.

None

Returns:

Name Type Description
dydt ndarray

Array of the time derivative of the state variables x, tau, and epsilon, concatenated into a 1D array.

References

[1] L. Chen, C. Li, T. Huang, H. G. Ahmad and Y. Chen, Physics Letters A, 2014, 378, 2924-2930

Source code in mnns/models.py
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
317
318
319
320
321
322
323
324
325
326
327
def SLT_HP_model(
    t: float,
    y: npt.NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: Optional[dict] = None,
) -> npt.NDArray:
    """
    SLT HP Model [1]. Provides the time derivative of the state variable `x`
    (the dimensionless version of `w`), `tau`, and `epsilon`. Assumes that
    these state variable are concatenated into a single 1D array input.

    Assumes voltage sources are used.

    Requires `NWN.graph["tau"]` to be set to the decay constant value.

    Parameters
    ----------
    t : float
        Current time to solve at.

    y : ndarray
        Array containing the state variables `x`, `tau`, and `epsilon` for each
        junction in the NWN concatenated into a single 1D array input.

    NWN : NanowireNetwork
        Input nanowire network graph.

    source_node : NWNNode or list of NWNNode
        Source node(s) of the input voltage.

    drain_node : NWNNode or list of NWNNode
        Drain/grounded node(s).

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

    window_func : Callable
        Function which inputs the state variable `x` as an array and returns
        the window function value as an array.

    solver : str
        SciPy sparse matrix equation solver.

    **kwargs
        Keyword arguments to pass to the SciPy sparse matrix equation solver.

    Returns
    -------
    dydt : ndarray
        Array of the time derivative of the state variables `x`, `tau`, and
        `epsilon`, concatenated into a 1D array.

    References
    ----------
    [1] L. Chen, C. Li, T. Huang, H. G. Ahmad and Y. Chen, *Physics Letters A*,
        2014, **378**, 2924-2930

    """
    if kwargs is None:
        kwargs = dict()

    # Unpack values
    w, tau, epsilon = np.split(y, 3)
    sigma, theta, a = NWN.graph["sigma"], NWN.graph["theta"], NWN.graph["a"]

    # Update all wire junction resistances
    R = NWN.update_resistance(w)

    # Find applied voltage at the current time
    applied_V = voltage_func(t)

    # Solve for voltage at each node
    *V, I = solve_network(
        NWN, source_node, drain_node, applied_V, "voltage", solver, **kwargs
    )
    V = np.asarray(V)

    # Get start and end indices
    start, end = NWN.wire_junction_indices()

    # Find voltage differences
    v0 = V[start]
    v1 = V[end]
    V_delta = np.abs(v0 - v1) * np.sign(applied_V)

    # Find derivatives
    l = V_delta / R
    dw_dt = (l - ((w - epsilon) / tau)) * window_func(w)
    dtau_dt = theta * l * (a - w)
    deps_dt = sigma * l * window_func(w)
    dydt = np.hstack((dw_dt, dtau_dt, deps_dt))

    return dydt

decay_HP_model

decay_HP_model(
    t: float,
    x: NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: dict = None,
) -> npt.NDArray

Decay HP Model [1]. Provides the time derivative of the state variable x (the dimensionless version of w). Assumes voltage sources are used.

Requires NWN.graph["tau"] to be set to the decay constant value.

Parameters:

Name Type Description Default
t float

Current time to solve at.

required
x ndarray

Array containing the state variable x for each junction in the NWN.

required
NWN NanowireNetwork

Input nanowire network graph.

required
source_node NWNNode or list of NWNNode

Source node(s) of the input voltage.

required
drain_node NWNNode or list of NWNNode

Drain/grounded node(s).

required
voltage_func Callable

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

required
window_func Callable

Function which inputs the state variable x as an array and returns the window function value as an array.

required
solver str

SciPy sparse matrix equation solver.

'spsolve'
**kwargs dict

Keyword arguments to pass to the SciPy sparse matrix equation solver.

None

Returns:

Name Type Description
dxdt ndarray

Array of the time derivative of the state variable x.

References

[1] H. O. Sillin, R. Aguilera, H.-H. Shieh, A. V. Avizienis, M. Aono, A. Z. Stieg and J. K. Gimzewski, Nanotechnology, 2013, 24, 384004.

Source code in mnns/models.py
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
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
def decay_HP_model(
    t: float,
    x: npt.NDArray,
    NWN: NanowireNetwork,
    source_node: NWNNode | list[NWNNode],
    drain_node: NWNNode | list[NWNNode],
    voltage_func: Callable,
    window_func: Callable,
    solver: str = "spsolve",
    kwargs: dict = None,
) -> npt.NDArray:
    """
    Decay HP Model [1]. Provides the time derivative of the state variable `x`
    (the dimensionless version of `w`). Assumes voltage sources are used.

    Requires `NWN.graph["tau"]` to be set to the decay constant value.

    Parameters
    ----------
    t : float
        Current time to solve at.

    x : ndarray
        Array containing the state variable `x` for each junction in the NWN.

    NWN : NanowireNetwork
        Input nanowire network graph.

    source_node : NWNNode or list of NWNNode
        Source node(s) of the input voltage.

    drain_node : NWNNode or list of NWNNode
        Drain/grounded node(s).

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

    window_func : Callable
        Function which inputs the state variable `x` as an array and returns
        the window function value as an array.

    solver : str
        SciPy sparse matrix equation solver.

    **kwargs
        Keyword arguments to pass to the SciPy sparse matrix equation solver.

    Returns
    -------
    dxdt : ndarray
        Array of the time derivative of the state variable `x`.

    References
    ----------
    [1] H. O. Sillin, R. Aguilera, H.-H. Shieh, A. V. Avizienis, M. Aono,
        A. Z. Stieg and J. K. Gimzewski, *Nanotechnology*, 2013, **24**, 384004.

    """
    if kwargs is None:
        kwargs = dict()

    # Update all wire junction resistances
    R = NWN.update_resistance(x)

    # Find applied voltage at the current time
    applied_V = voltage_func(t)

    # Solve for voltage at each node
    *V, I = solve_network(
        NWN, source_node, drain_node, applied_V, "voltage", solver, **kwargs
    )
    V = np.array(V)

    # Get start and end indices
    start, end = NWN.wire_junction_indices()

    # Find voltage differences
    v0 = V[start]
    v1 = V[end]
    V_delta = np.abs(v0 - v1) * np.sign(applied_V)

    # Get decay constant
    tau = NWN.graph["tau"]

    # Find dw/dt
    dxdt = (V_delta / R * window_func(x)) - (x / tau)

    return dxdt

resist_func

resist_func(
    NWN: NanowireNetwork, w: float | NDArray
) -> float | npt.NDArray

Linear resistance function in nondimensionalized form.

Obtained from Strukov et al. Nature, 2008, 453, 80-83.

Parameters:

Name Type Description Default
NWN NanowireNetwork

Nanowire network.

required
w ndarray or scalar

Nondimensionalized state variable of the memristor element(s).

required

Returns:

Name Type Description
R ndarray or scalar

Resistance of the memristor element(s).

Source code in mnns/models.py
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
def resist_func(
    NWN: NanowireNetwork, w: float | npt.NDArray
) -> float | npt.NDArray:
    """
    Linear resistance function in nondimensionalized form.

    Obtained from Strukov et al. *Nature*, 2008, **453**, 80-83.

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

    w : ndarray or scalar
        Nondimensionalized state variable of the memristor element(s).

    Returns
    -------
    R : ndarray or scalar
        Resistance of the memristor element(s).

    """
    Roff_Ron = NWN.graph["units"]["Roff_Ron"]
    R = w * (1 - Roff_Ron) + Roff_Ron
    return R