From 3547db353dfac9756e1efd4dd22e5d5b4208b627 Mon Sep 17 00:00:00 2001 From: LemurPwned Date: Thu, 30 Jan 2025 22:30:08 +0100 Subject: [PATCH] small updates to resistances and cvector + noise fixes --- cmtj/utils/resistance.py | 41 ++++++++++++++++++++++++++++++++++++++ core/cvector.hpp | 5 +++++ core/noise.hpp | 5 +++-- docs/api/core.md | 2 +- docs/index.md | 43 ++++++++++++++++++++++++++++++++++++++-- python/cmtj.cpp | 5 +++-- 6 files changed, 94 insertions(+), 7 deletions(-) diff --git a/cmtj/utils/resistance.py b/cmtj/utils/resistance.py index 10e624c..3c5b9bd 100644 --- a/cmtj/utils/resistance.py +++ b/cmtj/utils/resistance.py @@ -274,6 +274,47 @@ def Rxx_parallel_bilayer_expr(): return Rlin_func, R_func +def GMR_expr(): + """Get the symbolic expression for the GMR. + :returns: GMR function + """ + GMR_s = sym.Symbol(r"\mathrm{GMR}") + theta1 = sym.Symbol(r"\theta_1") + phi1 = sym.Symbol(r"\phi_1") + m1 = sym.Matrix( + [ + sym.sin(theta1) * sym.cos(phi1), + sym.sin(theta1) * sym.sin(phi1), + sym.cos(theta1), + ] + ) + theta2 = sym.Symbol(r"\theta_2") + phi2 = sym.Symbol(r"\phi_2") + m2 = sym.Matrix( + [ + sym.sin(theta2) * sym.cos(phi2), + sym.sin(theta2) * sym.sin(phi2), + sym.cos(theta2), + ] + ) + Rf = GMR_s * (1 - m1.dot(m2)) / 2 + dRdt1 = sym.diff(Rf, theta1) + dRdp1 = sym.diff(Rf, phi1) + dRdt2 = sym.diff(Rf, theta2) + dRdp2 = sym.diff(Rf, phi2) + linearised_terms = sym.symbols(r"\partial\theta_1, \partial\phi_1, \partial\theta_2, \partial\phi_2") + dRf = ( + dRdt1 * linearised_terms[0] + + dRdp1 * linearised_terms[1] + + dRdt2 * linearised_terms[2] + + dRdp2 * linearised_terms[3] + ) + + Rf_func = sym.lambdify([GMR_s, [theta1, phi1, theta2, phi2]], Rf) + dRf_func = sym.lambdify([GMR_s, [theta1, phi1, theta2, phi2], linearised_terms], dRf) + return Rf_func, dRf_func + + def Rxx_series_bilayer_expr(): """Get the symbolic expressions for the series and linearised resistance of a bilayer system. diff --git a/core/cvector.hpp b/core/cvector.hpp index 301a207..90271e1 100644 --- a/core/cvector.hpp +++ b/core/cvector.hpp @@ -212,6 +212,11 @@ template class CVector { ss << "[x:" << this->x << ", y:" << this->y << ", z:" << this->z << "]"; return ss.str(); } + + static CVector fromSpherical(T theta, T phi, T r = 1.0) { + return CVector(r * sin(theta) * cos(phi), r * sin(theta) * sin(phi), + r * cos(theta)); + } }; #endif // CORE_CVECTOR_HPP_ diff --git a/core/noise.hpp b/core/noise.hpp index 689f447..cba0397 100644 --- a/core/noise.hpp +++ b/core/noise.hpp @@ -134,7 +134,7 @@ template class OneFNoise { } }; -std::mt19937 generator(std::random_device{}()); +// std::mt19937 generator(std::random_device{}()); template class BufferedAlphaNoise : public NullTicker { protected: std::vector> bufferWhite, bufferColoured; @@ -148,6 +148,7 @@ template class BufferedAlphaNoise : public NullTicker { inv; // configs for forward and inverse real fft unsigned int internalCounter = 0; unsigned int refills = 0; + std::mt19937 generator; public: /** @@ -160,7 +161,7 @@ template class BufferedAlphaNoise : public NullTicker { */ BufferedAlphaNoise(unsigned int bufferSize, T alpha, T std, T scale) : bufferSize(bufferSize), alpha(alpha), scale(scale) { - + this->generator = std::mt19937(std::random_device{}()); this->bufferColoured.resize(2 * bufferSize); this->bufferWhite.resize(2 * bufferSize); this->result.resize(bufferSize); diff --git a/docs/api/core.md b/docs/api/core.md index f7992cf..8d42ec8 100644 --- a/docs/api/core.md +++ b/docs/api/core.md @@ -1 +1 @@ -::: cmtj +::: cmtj.core diff --git a/docs/index.md b/docs/index.md index c87648d..a073cd2 100644 --- a/docs/index.md +++ b/docs/index.md @@ -6,7 +6,7 @@ [![pages-build-deployment](https://github.com/LemurPwned/cmtj/actions/workflows/pages/pages-build-deployment/badge.svg?branch=gh-pages)](https://github.com/LemurPwned/cmtj/actions/workflows/pages/pages-build-deployment) [![Version](https://img.shields.io/pypi/v/cmtj)](https://pypi.org/project/cmtj/) [![License](https://img.shields.io/pypi/l/cmtj.svg)](https://github.com/LemurPwned/cmtj/blob/master/LICENSE) -[![Streamlit](https://static.streamlit.io/badges/streamlit_badge_black_white.svg)](http://cmtj-simulations.streamlit.app/) +[![Streamlit](https://static.streamlit.io/badges/streamlit_badge_black_white.svg)](https://cmtj-app.streamlit.app/spectrum) ![Downloads](https://img.shields.io/pypi/dm/cmtj.svg) ## Table of contents @@ -14,6 +14,7 @@ - [CMTJ](#cmtj) - [Table of contents](#table-of-contents) - [Short description](#short-description) + - [What can you simulate?](#what-can-you-simulate) - [Web GUI](#web-gui) - [Quickstart](#quickstart) - [Installation :rocket:](#installation-rocket) @@ -34,9 +35,47 @@ The `cmtj` name may be misleading -- the MTJ (Magnetic Tunnel Junctions) are not The library allows for macromagnetic simulation of various multilayer spintronic structures. The package uses C++ implementation of (s)LLGS (stochastic Landau-Lifschitz-Gilbert-Slonczewski) equation with various field contributions included for instance: anisotropy, interlayer exchange coupling, demagnetisation, dipole fields etc. It is also possible to connect devices in parallel or in series to have electrically coupled arrays. +### What can you simulate? + +Below is a brief list of examples (it's not exhaustive! Check the docs for more). + +**Magnetic devices:** + +- Magnetic Tunnel Junctions + - Voltage-Driven Magnetic Tunnel Junctions + - Spin-Torque Oscillators + - VCMA sensors and devices + - Magnetic Tunnel Junction Arrays +- SOT devices + - Current-Driven SOT +- Advanced device coupling +- Reservoirs (dipole coupling) +- Electrically coupled MTJs +- Base equations + - Landau-Lifshitz-Gilbert-Slonczewski equation + - Stochastic Landau-Lifshitz-Gilbert-Slonczewski equation + - Landau-Lifshitz-Gilbert-Bloch equation +- Domain wall motion + +**Experimental methods:** + +Some of the experimental methods available: + +- PIMM +- Spin-Diode +- CIMS +- R(H), M(H) + ## Web GUI -Check out the [streamlit hosted demo here](http://cmtj-simulations.streamlit.app/). You can simulate PIMM spectra and Spin-Diode spectra there. Let us know if you have any issues with the demo. +Check out the [streamlit hosted demo here](https://cmtj-app.streamlit.app/spectrum). +You can simulate: + +- PIMM spectra and Spin-Diode spectra +- Try some optimization fitting +- Fit multi-domain or multi-level M(H) or R(H) loops in [Domain mode](https://cmtj-app.streamlit.app) + +Let us know if you have any issues with the demo. ## Quickstart diff --git a/python/cmtj.cpp b/python/cmtj.cpp index 7b12417..6a68fa0 100644 --- a/python/cmtj.cpp +++ b/python/cmtj.cpp @@ -119,7 +119,8 @@ PYBIND11_MODULE(cmtj, m) { [](const DVector& v, const int key) { return v[key]; }) .def("__len__", [](const DVector& v) { return 3; }) .def("__str__", py::overload_cast<>(&DVector::toString)) - .def("__repr__", py::overload_cast<>(&DVector::toString)); + .def("__repr__", py::overload_cast<>(&DVector::toString)) + .def_static("fromSpherical", &DVector::fromSpherical, "theta"_a, "phi"_a, "r"_a = 1.0); py::implicitly_convertible, DVector>(); py::implicitly_convertible, DVector>(); @@ -240,7 +241,7 @@ PYBIND11_MODULE(cmtj, m) { .def_readonly("cellSurface", &DLayer::cellSurface) .def_readonly("demagTensor", &DLayer::demagTensor) // noise - .def("setAlphaNoise", &DLayer::setAlphaNoise) + .def("setAlphaNoise", &DLayer::setAlphaNoise, "alpha"_a, "std"_a, "scale"_a, "axis"_a = Axis::all) .def("setOneFNoise", &DLayer::setOneFNoise) // getters .def("getId", &DLayer::getId)