{
 "cells": [
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "b01c5c92-e134-42f9-bf8e-5c3f5b553ef7",
   "metadata": {},
   "outputs": [],
   "source": [
    "# pyright: reportUnknownArgumentType=false\n",
    "from rich.theme import Theme\n",
    "from rich.console import Console\n",
    "\n",
    "from finesse.model import Model\n",
    "from finesse.analysis.actions import (\n",
    "    TemporaryParameters,\n",
    "    Change,\n",
    "    Maximize,\n",
    "    Minimize,\n",
    "    Series,\n",
    "    FrequencyResponse,\n",
    "    Xaxis,\n",
    "    Noxaxis,\n",
    ")\n",
    "from finesse.solutions import SeriesSolution\n",
    "\n",
    "from matplotlib.pyplot import figure, show\n",
    "from matplotlib.axes import Axes\n",
    "\n",
    "from numpy import geomspace, linspace, sqrt\n",
    "from science_signal import Signal\n",
    "from scipy.io.matlab import loadmat\n",
    "\n",
    "from numpy.typing import NDArray\n",
    "from typing import Any, NamedTuple, Literal, Callable\n",
    "\n",
    "from pathlib import Path"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "f0b4199a-b9b0-4969-97c7-497faa662b94",
   "metadata": {},
   "outputs": [],
   "source": [
    "from gettext import install\n",
    "from logging import getLogger"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "71f985ab-ce47-4621-86c5-d4dc41a267d4",
   "metadata": {},
   "outputs": [],
   "source": [
    "install(__name__)\n",
    "logger = getLogger(__name__)\n",
    "theme = Theme(\n",
    "    {\n",
    "        \"strong\": \"cyan underline\",\n",
    "        \"result\": \"red bold\",\n",
    "        \"error\": \"red underline bold\",\n",
    "    }\n",
    ")\n",
    "console = Console(theme=theme)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "3b5dbd67-3baa-4cba-a5a2-368d306735a5",
   "metadata": {},
   "outputs": [],
   "source": [
    "C_DARK_FRINGE = 8e-3"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "d6999b95-7956-4d2a-a97e-93f3d12f5c74",
   "metadata": {},
   "outputs": [],
   "source": [
    "%matplotlib ipympl\n",
    "model_file = Path(\"model.kat\")\n",
    "model = Model()\n",
    "model.phase_config(zero_k00=False, zero_tem00_gouy=True)\n",
    "model.modes(modes=\"off\")  # pyright: ignore[reportUnusedCallResult]\n",
    "model.parse(model_file.read_text())\n",
    "model.lambda0 = model.get(\"wavelength\")\n",
    "model.plot_graph()  # pyright: ignore[reportUnusedCallResult]\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "d2b2ac68-0730-4f68-a4b1-ca036eadc880",
   "metadata": {},
   "outputs": [],
   "source": [
    "result = model.run(\n",
    "    TemporaryParameters(\n",
    "        Series(\n",
    "            Change(\n",
    "                {\n",
    "                    \"SR.misaligned\": True,\n",
    "                    \"PR.misaligned\": True,\n",
    "                    \"eom1.midx\": 0,\n",
    "                    \"eom2.midx\": 0,\n",
    "                    \"eom3.midx\": 0,\n",
    "                    \"eom4.midx\": 0,\n",
    "                }\n",
    "            ),\n",
    "            Maximize(\n",
    "                model.get(\"NE_p1\"),\n",
    "                model.get(\"NORTH_ARM.DC\"),\n",
    "                bounds=[-180, 180],\n",
    "                tol=1e-14,\n",
    "            ),\n",
    "            Maximize(\n",
    "                model.get(\"WE_p1\"),\n",
    "                model.get(\"WEST_ARM.DC\"),\n",
    "                bounds=[-180, 180],\n",
    "                tol=1e-14,\n",
    "            ),\n",
    "            Minimize(\n",
    "                model.get(\"SR_p2\"), model.get(\"MICH.DC\"), bounds=[-180, 180], tol=1e-14\n",
    "            ),\n",
    "            Change(\n",
    "                {\n",
    "                    \"PR.misaligned\": False,\n",
    "                }\n",
    "            ),\n",
    "            Maximize(\n",
    "                model.get(\"PR_p2\"), model.get(\"PRCL.DC\"), bounds=[-180, 180], tol=1e-14\n",
    "            ),\n",
    "            Change(\n",
    "                {\n",
    "                    \"SR.misaligned\": False,\n",
    "                }\n",
    "            ),\n",
    "            Maximize(\n",
    "                model.get(\"B1_DC\"), model.get(\"SRCL.DC\"), bounds=[-180, 180], tol=1e-14\n",
    "            ),\n",
    "            Change(\n",
    "                {\n",
    "                    \"SRCL.DC\": -90,\n",
    "                },\n",
    "                relative=True,\n",
    "            ),\n",
    "        ),\n",
    "        exclude=[\n",
    "            \"NE.phi\",\n",
    "            \"NI.phi\",\n",
    "            \"WE.phi\",\n",
    "            \"WI.phi\",\n",
    "            \"SR.phi\",\n",
    "            \"PR.phi\",\n",
    "            \"NORTH_ARM.DC\",\n",
    "            \"WEST_ARM.DC\",\n",
    "            \"DARM.DC\",\n",
    "            \"MICH.DC\",\n",
    "            \"PRCL.DC\",\n",
    "            \"SRCL.DC\",\n",
    "            \"SR.misaligned\",\n",
    "            \"eom1.midx\",\n",
    "            \"eom2.midx\",\n",
    "            \"eom3.midx\",\n",
    "            \"eom4.midx\",\n",
    "        ],\n",
    "    ),\n",
    ")\n",
    "model._settings.phase_config.zero_k00 = False\n",
    "model.fsig.f = 1"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "7504821e-b33d-4396-b86c-06e90477eb8e",
   "metadata": {},
   "outputs": [],
   "source": [
    "def compute_solutions(\n",
    "    model: Model, DOF: str, padding: float, nb: int = 10000\n",
    ") -> SeriesSolution:\n",
    "    return model.run(\n",
    "        Xaxis(\n",
    "            model.get(DOF).DC,\n",
    "            \"lin\",\n",
    "            model.get(DOF).DC - padding,\n",
    "            model.get(DOF).DC + padding,\n",
    "            nb,\n",
    "        )\n",
    "    )\n",
    "\n",
    "\n",
    "def display_ax(\n",
    "    ax: Axes,\n",
    "    solution: SeriesSolution,\n",
    "    model: Model,\n",
    "    DOF: str,\n",
    "    padding: float,\n",
    "    nb: int = 10000,\n",
    ") -> Axes:\n",
    "    x = linspace(model.get(DOF).DC - padding, model.get(DOF).DC + padding, nb + 1)\n",
    "    _ = ax.semilogy(x, solution[\"SR_p2\"], label=\"dark fringe\")\n",
    "    _ = ax.semilogy(x, solution[\"NE_p1\"], label=\"north cavity\")\n",
    "    _ = ax.semilogy(x, solution[\"WE_p1\"], label=\"west cavity\")\n",
    "    _ = ax.vlines(\n",
    "        [model.get(DOF).DC],\n",
    "        min(solution[\"SR_p2\"]),\n",
    "        max(solution[\"NE_p1\"]),\n",
    "        colors=\"red\",\n",
    "    )\n",
    "    _ = ax.set_ylabel(\"power (W)\")\n",
    "    ax.grid()\n",
    "    _ = ax.legend()\n",
    "    return ax\n",
    "\n",
    "\n",
    "class DisplayData(NamedTuple):\n",
    "    DOF: str\n",
    "    padding: float\n",
    "\n",
    "\n",
    "data: list[DisplayData] = [\n",
    "    DisplayData(\"NORTH_ARM\", 10),\n",
    "    DisplayData(\"WEST_ARM\", 10),\n",
    "    DisplayData(\"PRCL\", 10),\n",
    "    DisplayData(\"MICH\", 10),\n",
    "    DisplayData(\"DARM\", 10),\n",
    "    DisplayData(\"CARM\", 10),\n",
    "]\n",
    "\n",
    "Figure = figure(figsize=(13, 10))\n",
    "nb = int(1e4)\n",
    "\n",
    "for i in range(len(data)):\n",
    "    element: DisplayData = data[i]\n",
    "    ax = Figure.add_subplot(3, 2, i + 1)\n",
    "    solution = compute_solutions(model, element.DOF, element.padding, nb)\n",
    "    _ = display_ax(ax, solution, model, element.DOF, element.padding, nb).set_xlabel(\n",
    "        \"{} value\".format(element.DOF)\n",
    "    )\n",
    "show()\n",
    "\n",
    "solution = model.run(Noxaxis())\n",
    "result = solution[\"B1_DC\"]\n",
    "start, stop, nb = 0, 1, 0\n",
    "while (abs(result - C_DARK_FRINGE) > 1e-4) and (nb < 100):\n",
    "    nb += 1\n",
    "    temp = start + (stop - start) / 2\n",
    "\n",
    "    model.DARM.DC = temp\n",
    "    solution = model.run(Noxaxis())\n",
    "    result = solution[\"B1_DC\"]\n",
    "    if result > C_DARK_FRINGE:\n",
    "        stop = temp\n",
    "    else:\n",
    "        start = temp\n",
    "console.print(\n",
    "    \"Degré de liberté [result]{dof}[/result] trouvé en [strong]{nb} pas[/strong] pour avoir une puissance de [result]{result} W[/result] sur B1\".format(\n",
    "        nb=nb, dof=model.DARM.DC, result=result\n",
    "    )\n",
    ")"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "3fbf07f1-22bc-445b-95ad-18a30adcc85d",
   "metadata": {},
   "outputs": [],
   "source": [
    "def show_evolution(\n",
    "    parameter: str,\n",
    "    model: Model,\n",
    "    values: NDArray[Any],\n",
    "    TFs: list[str],\n",
    "    phase: int | float = 45,\n",
    "):\n",
    "    power_detector = \"B1.I\"\n",
    "    model = model.deepcopy()\n",
    "    model.SNEB.phi = model.NE.phi - phase\n",
    "    model.SWEB.phi = model.WE.phi - phase\n",
    "    model.SDB1.phi = model.SR.phi + phase\n",
    "\n",
    "    if len(TFs) == 0:\n",
    "        console.print(\"[error]Nothing to show[/error]\")\n",
    "\n",
    "    Figure = figure(figsize=(7, 5 * len(TFs)))\n",
    "    Figure.suptitle(\"TF in function of {}\".format(parameter))\n",
    "\n",
    "    for i in range(len(TFs)):\n",
    "        _ = Figure.add_subplot(len(TFs), 1, i + 1)\n",
    "\n",
    "    temp_value = model.get(parameter).eval()\n",
    "    for value in values:\n",
    "        index = 0\n",
    "        model.set(parameter, value)\n",
    "\n",
    "        DARM = model.run(\n",
    "            FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [power_detector])\n",
    "        )\n",
    "        for bench in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
    "            if bench in TFs:\n",
    "                result = model.run(\n",
    "                    FrequencyResponse(\n",
    "                        geomspace(5, 10000, 1000),\n",
    "                        [\"{}_z\".format(bench)],\n",
    "                        [power_detector],\n",
    "                    )\n",
    "                )\n",
    "                _ = Figure.get_axes()[index].set_title(bench)\n",
    "                _ = Figure.get_axes()[index].loglog(\n",
    "                    result.f,\n",
    "                    abs(result[power_detector, \"{}_z\".format(bench)])\n",
    "                    / abs(DARM[power_detector, \"DARM\"])\n",
    "                    / model.space_NI_NE.L.eval(),\n",
    "                    label=\"{} = {:.2E}\".format(parameter, value),\n",
    "                )\n",
    "                _ = Figure.get_axes()[index].set_xlabel(\"Frequencies (Hz)\")\n",
    "                _ = Figure.get_axes()[index].set_ylabel(\"$\\\\frac{ m } { m }$\")\n",
    "                if phase == 45:\n",
    "                    _ = Figure.get_axes()[index].set_title(\"{} $K_P$\".format(bench))\n",
    "                else:\n",
    "                    _ = Figure.get_axes()[index].set_title(\"{} $K_n$\".format(bench))\n",
    "\n",
    "                index += 1\n",
    "    for i in range(len(TFs)):\n",
    "        _ = Figure.get_axes()[i].grid(True, \"both\", \"both\")\n",
    "        _ = Figure.get_axes()[i].legend()\n",
    "    show()\n",
    "    model.set(parameter, temp_value)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "b2b8bc8d-dd73-4a45-af9d-179b3f7a3c7f",
   "metadata": {},
   "outputs": [],
   "source": [
    "show_evolution(\"NE.T\", model, geomspace(1e-8, 1.5e-5, 10), [\"SNEB\"])\n",
    "show_evolution(\"NE.T\", model, geomspace(1e-8, 1.5e-5, 10), [\"SNEB\"], 0)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "4a514ded-0777-41a7-82ce-380894d7be44",
   "metadata": {},
   "outputs": [],
   "source": [
    "show_evolution(\"WE.T\", model, geomspace(1e-8, 1.5e-5, 10), [\"SWEB\"])\n",
    "show_evolution(\"WE.T\", model, geomspace(1e-8, 1.5e-5, 10), [\"SWEB\"], 0)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "5a7d30d4-4d9f-4715-a1d6-a0330e8450a7",
   "metadata": {},
   "outputs": [],
   "source": [
    "# show_evolution(\"SR.T\", model, linspace(0.30, 0.50, 10), [\"SDB1\"])"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "c7dca854-5dff-483f-8a72-efaf0b80d87d",
   "metadata": {},
   "outputs": [],
   "source": [
    "# show_evolution(\"NI.T\", model, linspace(1.35e-2, 1.39e-2, 10), [\"SNEB\", \"SWEB\", \"SDB1\"])"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "b31a0f43-0a07-44b4-9264-b5001350eb9d",
   "metadata": {},
   "outputs": [],
   "source": [
    "# show_evolution(\"WI.T\", model, linspace(1.35e-2, 1.39e-2, 10), [\"SNEB\", \"SWEB\", \"SDB1\"])"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "fbb27d21-7f79-4021-8833-98a96a76ce18",
   "metadata": {},
   "outputs": [],
   "source": [
    "def compare_TF(\n",
    "    model: Model,\n",
    "    bench: str,\n",
    "    phase: Literal[45] | Literal[0],\n",
    ") -> tuple[Signal, Signal, Signal]:\n",
    "    power_detector = \"B1.I\"\n",
    "    model = model.deepcopy()\n",
    "    model.SNEB.phi = model.NE.phi - phase\n",
    "    model.SWEB.phi = model.WE.phi - phase\n",
    "    model.SDB1.phi = model.SR.phi + phase\n",
    "\n",
    "    modelisation_data = loadmat(Path(\"optickle.mat\"))\n",
    "\n",
    "    DARM = model.run(\n",
    "        FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [power_detector]),\n",
    "    )\n",
    "    TF = model.run(\n",
    "        FrequencyResponse(\n",
    "            geomspace(5, 10000, 1000), [\"{}_z\".format(bench)], [power_detector]\n",
    "        ),\n",
    "    )\n",
    "\n",
    "    TF_finesse = Signal(\n",
    "        TF.f,\n",
    "        abs(TF[power_detector, \"{}_z\".format(bench)])\n",
    "        / abs(DARM[power_detector, \"DARM\"])\n",
    "        / model.space_NI_NE.L,\n",
    "    )\n",
    "    TF_optickle = Signal(\n",
    "        modelisation_data[\"freq\"][0],\n",
    "        abs(modelisation_data[\"{}coupling\".format(bench)][int(phase / 45)]),  # 1 or 0\n",
    "    )\n",
    "    return abs(TF_finesse - TF_optickle)/(abs(TF_finesse) + abs(TF_optickle)), TF_finesse, TF_optickle"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "81659108-bcd9-4465-b733-b7cf82719f26",
   "metadata": {},
   "outputs": [],
   "source": [
    "def compare_allTF(\n",
    "    model: Model,\n",
    "    bench: str,\n",
    ") -> tuple[Signal, Signal, Signal]:\n",
    "    power_detector = \"B1.I\"\n",
    "    model = model.deepcopy()\n",
    "    model.SNEB.phi = model.NE.phi - 0\n",
    "    model.SWEB.phi = model.WE.phi - 0\n",
    "    model.SDB1.phi = model.SR.phi + 0\n",
    "\n",
    "    modelisation_data = loadmat(Path(\"optickle.mat\"))\n",
    "\n",
    "    in_DARM = model.run(\n",
    "        FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [power_detector]),\n",
    "    )\n",
    "    in_TF = model.run(\n",
    "        FrequencyResponse(\n",
    "            geomspace(5, 10000, 1000), [\"{}_z\".format(bench)], [power_detector]\n",
    "        ),\n",
    "    )\n",
    "    model.SNEB.phi = model.NE.phi - 45\n",
    "    model.SWEB.phi = model.WE.phi - 45\n",
    "    model.SDB1.phi = model.SR.phi + 45\n",
    "\n",
    "    quad_DARM = model.run(\n",
    "        FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [power_detector]),\n",
    "    )\n",
    "    quad_TF = model.run(\n",
    "        FrequencyResponse(\n",
    "            geomspace(5, 10000, 1000), [\"{}_z\".format(bench)], [power_detector]\n",
    "        ),\n",
    "    )\n",
    "\n",
    "    TF_finesse = Signal(\n",
    "        in_TF.f,\n",
    "        sqrt((abs(in_TF[power_detector, \"{}_z\".format(bench)])\n",
    "        / abs(in_DARM[power_detector, \"DARM\"])\n",
    "        / model.space_NI_NE.L)**2 + (abs(quad_TF[power_detector, \"{}_z\".format(bench)])\n",
    "        / abs(quad_DARM[power_detector, \"DARM\"])\n",
    "        / model.space_NI_NE.L)**2),\n",
    "    )\n",
    "    TF_optickle = Signal(\n",
    "        modelisation_data[\"freq\"][0],\n",
    "        sqrt(abs(modelisation_data[\"{}coupling\".format(bench)][0]**2 + modelisation_data[\"{}coupling\".format(bench)][1]**2)),\n",
    "    )\n",
    "    return abs(TF_finesse - TF_optickle)/(abs(TF_finesse) + abs(TF_optickle)), TF_finesse, TF_optickle"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "172943ef-c784-4481-9d97-c7330c2a5bec",
   "metadata": {},
   "outputs": [],
   "source": [
    "Figure = figure(figsize = (14, 5))\n",
    "result, TF_finesse, TF_optickle = compare_TF(model, \"SNEB\", 45)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label = \"difference\")\n",
    "_ = Figure.gca().loglog(TF_finesse.x, TF_finesse.y, label = \"finesse\")\n",
    "_ = Figure.gca().loglog(TF_optickle.x, TF_optickle.y, label = \"optickle\")\n",
    "_ = Figure.gca().legend()\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "7669861e-c40d-4cdb-bbd9-e256d023ac6c",
   "metadata": {},
   "outputs": [],
   "source": [
    "Figure = figure(figsize=(7, 5))\n",
    "result, _, _ = compare_TF(model, \"SNEB\", 45)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SNEB\")\n",
    "result, _, _ = compare_TF(model, \"SWEB\", 45)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SWEB\")\n",
    "result, _, _ = compare_TF(model, \"SDB1\", 45)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SDB1\")\n",
    "_ = Figure.gca().legend()\n",
    "_ = Figure.gca().set_title(\"Difference between Optickle and Finesse Transfer Function's module ($K_P$)\")\n",
    "_ = Figure.gca().set_xlabel(\"Frequencies (Hz)\")\n",
    "_ = Figure.gca().set_ylabel(\"$\\\\frac { m } { m }$\")\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()\n",
    "Figure = figure(figsize=(7, 5))\n",
    "result, _, _ = compare_TF(model, \"SNEB\", 0)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SNEB\")\n",
    "result, _, _ = compare_TF(model, \"SWEB\", 0)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SWEB\")\n",
    "result, _, _ = compare_TF(model, \"SDB1\", 0)\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SDB1\")\n",
    "_ = Figure.gca().legend()\n",
    "_ = Figure.gca().set_title(\"Difference between Optickle and Finesse Transfer Function's module ($K_n$)\")\n",
    "_ = Figure.gca().set_xlabel(\"Frequencies (Hz)\")\n",
    "_ = Figure.gca().set_ylabel(\"$\\\\frac { m } { m }$\")\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "aaf61dc8-c3e0-4476-9bcb-ba85e8cc6b3c",
   "metadata": {},
   "outputs": [],
   "source": [
    "Figure = figure(figsize=(7, 5))\n",
    "result, _, _ = compare_allTF(model, \"SNEB\")\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SNEB\")\n",
    "result, _, _ = compare_allTF(model, \"SWEB\")\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SWEB\")\n",
    "result, _, _ = compare_allTF(model, \"SDB1\")\n",
    "_ = Figure.gca().loglog(result.x, result.y, label=\"SDB1\")\n",
    "_ = Figure.gca().legend()\n",
    "_ = Figure.gca().set_title(\"Difference between Optickle and Finesse Transfer Function's module\")\n",
    "_ = Figure.gca().set_xlabel(\"Frequencies (Hz)\")\n",
    "_ = Figure.gca().set_ylabel(\"$\\\\frac { m } { m }$\")\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "d43f1976-7694-418d-bf65-fe2921a500fa",
   "metadata": {},
   "outputs": [],
   "source": [
    "def change_parameter(\n",
    "    parameter: str,\n",
    "    model: Model,\n",
    "    values: NDArray[Any],\n",
    "    func: Callable[Any, Signal],\n",
    "    args: list[Any],\n",
    ") -> list[Signal]:\n",
    "    list_result = []\n",
    "    model = model.deepcopy()\n",
    "    for value in values:\n",
    "        model.set(parameter, value)\n",
    "        list_result.append(func(model, *args))\n",
    "    return list_result"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "6f6e3b88-e98a-4a14-8e71-b00ac28cca4b",
   "metadata": {},
   "outputs": [],
   "source": [
    "x = linspace(2e-6, 6e-6, 10)\n",
    "SDB1_quad = change_parameter(\"NE.T\", model, x, compare_TF, [\"SNEB\", 0])\n",
    "console.print(model.NI.T)\n",
    "Figure = figure(figsize=(14, 5))\n",
    "for i in range(len(x)):\n",
    "    _ = Figure.gca().loglog(SDB1_quad[i][0].x, SDB1_quad[i][0].y, label = \"{}\".format(x[i]))\n",
    "_ = Figure.gca().legend()\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "e4227a4f-6c82-41d8-b586-adcd634cdbd8",
   "metadata": {},
   "outputs": [],
   "source": [
    "x = linspace(4.29e-6, 4.31e-6, 10)\n",
    "SDB1_quad = change_parameter(\"WE.T\", model, x, compare_TF, [\"SWEB\", 0])\n",
    "console.print(model.WE.T)\n",
    "Figure = figure(figsize=(14, 5))\n",
    "for i in range(len(x)):\n",
    "    _ = Figure.gca().loglog(SDB1_quad[i][0].x, SDB1_quad[i][0].y, label = \"{}\".format(x[i]))\n",
    "_ = Figure.gca().legend()\n",
    "Figure.gca().grid(True, \"both\", \"both\")\n",
    "show()"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "26e77672-0cc6-4239-8d4e-985ce3396e7e",
   "metadata": {},
   "outputs": [],
   "source": [
    "console.log(model.NE.T)\n",
    "console.log(model.WE.T)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "bdacc0d4-2b93-4c40-b641-9d563bdb431a",
   "metadata": {},
   "outputs": [],
   "source": []
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3 (ipykernel)",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.13.2"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}