finesse-simulation-O4/main.ipynb
2025-04-01 14:47:27 +02:00

754 lines
26 KiB
Text

{
"cells": [
{
"cell_type": "markdown",
"id": "a6ba3eb0-8f27-4ebd-b407-3f25f449c6bf",
"metadata": {},
"source": [
"# Imports"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "bd9299aa-a531-468c-b04b-798b06315f41",
"metadata": {},
"outputs": [],
"source": [
"# pyright: reportUnknownArgumentType=false, reportCallIssue=false, reportAttributeAccessIssue=false, reportOptionalSubscript=false, reportArgumentType=false\n",
"from rich.console import Console\n",
"from rich.table import Table\n",
"from rich.theme import Theme\n",
"\n",
"from finesse.model import Model\n",
"from finesse.analysis.actions.axes import Noxaxis, Xaxis\n",
"from finesse.solutions import SeriesSolution\n",
"from finesse.analysis.actions import (\n",
" TemporaryParameters,\n",
" Change,\n",
" Maximize,\n",
" Minimize,\n",
" Series,\n",
" FrequencyResponse,\n",
")\n",
"from finesse.components import Mirror, SignalGenerator\n",
"from finesse.detectors import QuantumNoiseDetector\n",
"\n",
"from pathlib import Path\n",
"from typing import NamedTuple\n",
"import re\n",
"\n",
"from matplotlib.axes import Axes\n",
"from matplotlib.pyplot import figure, show\n",
"\n",
"\n",
"from numpy import linspace, geomspace, pi, angle, where, diff, mean, loadtxt, load\n",
"from scipy.io.matlab import loadmat\n",
"\n",
"from science_signal import Signal"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "4c038d40-1d01-49cb-9182-a9a0e94d0d40",
"metadata": {},
"outputs": [],
"source": [
"from gettext import install\n",
"from logging import getLogger"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "5d4f2612-c5ea-4b21-a326-7074022966bc",
"metadata": {},
"outputs": [],
"source": [
"install(__name__)\n",
"logger = getLogger(__name__)\n",
"theme = Theme(\n",
" {\n",
" \"strong\": \"cyan underline\",\n",
" \"result\": \"red bold\",\n",
" }\n",
")\n",
"console = Console(theme=theme)"
]
},
{
"cell_type": "markdown",
"id": "eb7d2340-c817-4309-9599-6d58070ff4ab",
"metadata": {},
"source": [
"## Paramètres généraux"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "8fc23eea-145e-4641-93e9-6f8989edca96",
"metadata": {},
"outputs": [],
"source": [
"C_POWER = 25 # en Whatt\n",
"C_DARK_FRINGE = 8e-3 # en Whatt"
]
},
{
"cell_type": "markdown",
"id": "3052aa2b-350e-4eb3-b31c-a4204ab84dac",
"metadata": {},
"source": [
"## Modèle simplifié de Virgo"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "d32480d0-8525-478a-9af5-b1a1c8b30f1d",
"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.laser.P = C_POWER\n",
"model.plot_graph() # pyright: ignore[reportUnusedCallResult]\n",
"show()"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "b742cd14-2149-437b-9194-249cf40849b1",
"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",
")"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "fbc0d68f-9b01-4d52-866c-18d49aedff32",
"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": "66c9be03-cedc-4d07-acbd-3269f3143cd4",
"metadata": {},
"outputs": [],
"source": [
"solution = model.run(Noxaxis())\n",
"console = Console()\n",
"table = Table(title=\"Puissances dans l'interferomètre\")\n",
"table.add_column(\"position\", justify=\"left\", style=\"white\")\n",
"table.add_column(\"puissance (W)\", justify=\"left\", style=\"cyan\")\n",
"\n",
"table.add_row(\"Injection\", str(model.get(\"laser\").P.eval()))\n",
"table.add_row(\"PR\", str(solution[\"PR_p1\"]))\n",
"table.add_row(\"cavité de recyclage de puissance\", str(solution[\"PR_p2\"]))\n",
"table.add_row(\"cavité ouest\", str(solution[\"WE_p1\"]))\n",
"table.add_row(\"cavité nord\", str(solution[\"NE_p1\"]))\n",
"table.add_row(\"frange noire\", str(solution[\"SR_p2\"]))\n",
"table.add_row(\"SNEB\", str(solution[\"SNEB_DC\"]))\n",
"table.add_row(\"SWEB\", str(solution[\"SWEB_DC\"]))\n",
"table.add_row(\"SDB1\", str(solution[\"SDB1_DC\"]))\n",
"\n",
"console.print(table)\n",
"\n",
"table = Table(title=\"DOF dans l'interferomètre\")\n",
"table.add_column(\"nom\", justify=\"left\", style=\"white\")\n",
"table.add_column(\"valeur\", justify=\"left\", style=\"magenta\")\n",
"\n",
"table.add_row(\"Bras nord\", str(model.get(\"NORTH_ARM.DC\")))\n",
"table.add_row(\"Bras ouest\", str(model.get(\"WEST_ARM.DC\")))\n",
"table.add_row(\"PR\", str(model.get(\"PRCL.DC\")))\n",
"table.add_row(\"SR\", str(model.get(\"SRCL.DC\")))\n",
"table.add_row(\"MICH\", str(model.get(\"MICH.DC\")))\n",
"\n",
"console.print(table)\n",
"\n",
"console = Console(theme=theme)\n",
"table = Table(title=\"\")\n",
"table.add_column(\"nom\", justify=\"left\", style=\"white\")\n",
"table.add_column(\"valeur\", justify=\"left\", style=\"cyan\")\n",
"for i in range(1, model.west_arm.info_parameter_table().table.shape[0]):\n",
" table.add_row(\n",
" str(model.west_arm.info_parameter_table().table[i, 0]),\n",
" str(model.west_arm.info_parameter_table().table[i, 1]),\n",
" )\n",
"console.print(table)"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "7fde067a-52b7-4798-8bd6-0890e27f33e2",
"metadata": {},
"outputs": [],
"source": [
"def get_QNLS(\n",
" model: Model, start: int = 5, stop: int = 1000, nb: int = 100\n",
") -> SeriesSolution:\n",
" new_model = model.deepcopy()\n",
" new_model.fsig.f = 1\n",
" new_model.add(SignalGenerator(\"darmx\", new_model.space_NI_NE.h, 1, 0))\n",
" new_model.add(SignalGenerator(\"darmy\", new_model.space_WI_WE.h, 1, 180))\n",
" new_model.add(QuantumNoiseDetector(\"NSR_with_RP\", new_model.SR.p2.o, True))\n",
" return new_model.run(Xaxis(new_model.get(\"fsig.f\"), \"log\", start, stop, nb))\n",
"\n",
"\n",
"model._settings.phase_config.zero_k00 = False\n",
"\n",
"solution = get_QNLS(model, 5, 5000, 10000)\n",
"\n",
"\n",
"def dumb_parse(value: str = \"\"):\n",
" regex = re.compile(\"\\\\((\\\\d+\\\\.\\\\d+e[+-]\\\\d{2})([+-]\\\\d+\\\\.\\\\d+e[+-]\\\\d{2})j\\\\)\")\n",
" result = re.search(regex, value)\n",
" if result:\n",
" return float(result.groups()[0]) + 1j * float(result.groups()[1])\n",
" raise Exception(value)\n",
"\n",
"\n",
"QNLS = load(\"sensitivities/finesse-virgo.npy\")\n",
"current_O4_sensitivity_ASD = loadtxt(\"sensitivities/O4_nominal_reference.txt\")\n",
"\n",
"Figure = figure(figsize=(14, 5))\n",
"_ = Figure.gca().loglog(\n",
" solution.x1, abs(solution[\"NSR_with_RP\"]), label=\"this lock process\"\n",
")\n",
"_ = Figure.gca().loglog(\n",
" QNLS[0],\n",
" QNLS[1],\n",
" label=\"packaged lock process\",\n",
")\n",
"_ = Figure.gca().loglog(\n",
" current_O4_sensitivity_ASD[0],\n",
" abs(current_O4_sensitivity_ASD[1]),\n",
" label=\"current nominal sensitivity during O4\",\n",
")\n",
"_ = Figure.gca().legend()\n",
"Figure.gca().grid(True, \"both\", \"both\")\n",
"show()\n",
"\n",
"solution = model.run(FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [\"B1.I\"]))\n",
"maximum_amplitude_step: float = max(abs(diff(angle(solution[\"B1.I\", \"DARM\"]))))\n",
"pole_index = round(\n",
" mean(\n",
" where(\n",
" abs(angle(solution[\"B1.I\", \"DARM\"]) + pi / 4) < maximum_amplitude_step * 2\n",
" )\n",
" )\n",
") # find the index where the curve is the closest to -45°\n",
"console.print(\n",
" \"Le [strong]pôle[/strong] de la fonction de transfert [strong]DARM[/strong] est à [result]{:.1f}[/result] Hz\".format(\n",
" solution.f[pole_index]\n",
" )\n",
")\n",
"\n",
"table = Table(title=\"Position des différents miroirs\")\n",
"table.add_column(\"miroir\", justify=\"left\", style=\"white\")\n",
"table.add_column(\"offset (°)\", justify=\"left\", style=\"white\")\n",
"table.add_column(\"offset (m)\", justify=\"left\", style=\"white\")\n",
"\n",
"for name in [\n",
" \"NE\",\n",
" \"NE_AR\",\n",
" \"NI\",\n",
" \"NI_AR\",\n",
" \"WE\",\n",
" \"WE_AR\",\n",
" \"WI\",\n",
" \"WI_AR\",\n",
" \"PR\",\n",
" \"PR_AR\",\n",
" \"SR\",\n",
" \"SR_AR\",\n",
"]:\n",
" element: Mirror = model.get(name)\n",
" table.add_row(\n",
" str(element.name),\n",
" str(element.phi.eval()),\n",
" str(element.phi.eval() * model.lambda0 / 180),\n",
" )\n",
"\n",
"console.print(table)"
]
},
{
"cell_type": "markdown",
"id": "fd5ff122-8d97-43f2-982a-084ee984b827",
"metadata": {},
"source": [
"## Comparaison avec Optickle"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "fd3b1078-0fd6-4e1b-9e55-b231d9464bdf",
"metadata": {},
"outputs": [],
"source": [
"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",
"B1_detector = \"B1.I\"\n",
"\n",
"quad_tf: dict[str, SeriesSolution] = dict()\n",
"in_tf: dict[str, SeriesSolution] = dict()\n",
"\n",
"for bench_name in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
" quad_tf[bench_name] = model.run(\n",
" FrequencyResponse(\n",
" geomspace(5, 10000, 1000), [\"{}_z\".format(bench_name)], [B1_detector]\n",
" )\n",
" )\n",
"\n",
"quad_tf[\"DARM\"] = model.run(\n",
" FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [B1_detector])\n",
")\n",
"\n",
"model.SNEB.phi = model.NE.phi\n",
"model.SWEB.phi = model.WE.phi\n",
"model.SDB1.phi = model.SR.phi\n",
"\n",
"for bench_name in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
" in_tf[bench_name] = model.run(\n",
" FrequencyResponse(\n",
" geomspace(5, 10000, 1000), [\"{}_z\".format(bench_name)], [B1_detector]\n",
" )\n",
" )\n",
"\n",
"in_tf[\"DARM\"] = model.run(\n",
" FrequencyResponse(geomspace(5, 10000, 1000), [\"DARM\"], [B1_detector])\n",
")\n",
"\n",
"modelisation_file = Path(\"optickle.mat\")\n",
"\n",
"modelisation_data = loadmat(modelisation_file)\n",
"coupling_data: dict[str, list[Signal]] = dict()\n",
"old_couplings = [\"SNEB\", \"SWEB\", \"SDB1\"]\n",
"for coupling in old_couplings:\n",
" coupling_data[coupling] = [\n",
" Signal(\n",
" modelisation_data[\"freq\"][0],\n",
" abs(values),\n",
" )\n",
" for values in modelisation_data[\"{}coupling\".format(coupling)]\n",
" ]\n",
"DARMcoupling = Signal(\n",
" modelisation_data[\"freq\"][0],\n",
" modelisation_data[\"DARMmat\"][0],\n",
")"
]
},
{
"cell_type": "markdown",
"id": "be43f8b2-eddf-4dfc-a342-a54017b571f6",
"metadata": {},
"source": [
"### En fonction de la phase"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "61c3d4e0-b8bc-48e0-83cd-5675c62feba1",
"metadata": {},
"outputs": [],
"source": [
"Figure = figure(figsize=(14, 5))\n",
"_ = Figure.suptitle(\"Comparaison du module des fonctions de transfert pour DARM\")\n",
"ax = Figure.add_subplot(1, 2, 1)\n",
"_ = ax.loglog(\n",
" quad_tf[\"DARM\"].f, abs(quad_tf[\"DARM\"][B1_detector, \"DARM\"]), label=\"Finesse\"\n",
")\n",
"_ = ax.loglog(DARMcoupling.x, abs(DARMcoupling.y), label=\"Optickle\")\n",
"_ = ax.set_title(\"En quadrature de phase\")\n",
"_ = ax.legend()\n",
"ax.grid(True, \"both\", \"both\")\n",
"ax = Figure.add_subplot(1, 2, 2)\n",
"_ = ax.loglog(in_tf[\"DARM\"].f, abs(in_tf[\"DARM\"][B1_detector, \"DARM\"]), label=\"Finesse\")\n",
"_ = ax.loglog(DARMcoupling.x, abs(DARMcoupling.y), label=\"Optickle\")\n",
"_ = ax.set_title(\"En phase\")\n",
"_ = ax.legend()\n",
"ax.grid(True, \"both\", \"both\")\n",
"\n",
"Figure = figure(figsize=(14, 5))\n",
"_ = Figure.suptitle(\"Comparaison de la phase des fonctions de transfert pour DARM\")\n",
"ax = Figure.add_subplot(1, 2, 1)\n",
"_ = ax.semilogx(\n",
" quad_tf[\"DARM\"].f,\n",
" angle(quad_tf[\"DARM\"][B1_detector, \"DARM\"]) * 180 / pi,\n",
" label=\"Finesse\",\n",
")\n",
"_ = ax.semilogx(\n",
" DARMcoupling.x, angle(DARMcoupling.y) * 180 / pi - 180, label=\"Optickle\"\n",
")\n",
"_ = ax.set_title(\"En quadrature de phase\")\n",
"_ = ax.hlines([-45], min(quad_tf[\"DARM\"].f), max(quad_tf[\"DARM\"].f), colors=\"red\")\n",
"_ = ax.legend()\n",
"ax.grid(True, \"both\", \"both\")\n",
"ax = Figure.add_subplot(1, 2, 2)\n",
"_ = ax.semilogx(\n",
" in_tf[\"DARM\"].f,\n",
" angle(in_tf[\"DARM\"][B1_detector, \"DARM\"]) * 180 / pi,\n",
" label=\"Finesse\",\n",
")\n",
"_ = ax.semilogx(\n",
" DARMcoupling.x, angle(DARMcoupling.y) * 180 / pi - 180, label=\"Optickle\"\n",
")\n",
"_ = ax.set_title(\"En phase\")\n",
"_ = ax.hlines([-45], min(quad_tf[\"DARM\"].f), max(quad_tf[\"DARM\"].f), colors=\"red\")\n",
"_ = ax.legend()\n",
"ax.grid(True, \"both\", \"both\")\n",
"\n",
"for bench_name in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
" in_index = 0\n",
" if bench_name == \"SDB1\":\n",
" in_index = 0\n",
" quad_index = (1 + in_index) % 2\n",
" Figure = figure(figsize=(14, 5))\n",
" _ = Figure.suptitle(\n",
" \"Comparaison du module des fonctions de transfert pour {}\".format(bench_name)\n",
" )\n",
" ax = Figure.add_subplot(1, 2, 1)\n",
" _ = ax.loglog(\n",
" quad_tf[bench_name].f,\n",
" abs(quad_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(quad_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"Finesse\",\n",
" )\n",
" _ = ax.loglog(\n",
" coupling_data[bench_name][quad_index].x,\n",
" abs(coupling_data[bench_name][quad_index].y),\n",
" label=\"Optickle\",\n",
" )\n",
" _ = ax.set_title(\"En Quadrature de phase\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" ax = Figure.add_subplot(1, 2, 2)\n",
" _ = ax.loglog(\n",
" in_tf[bench_name].f,\n",
" abs(in_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(in_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"Finesse\",\n",
" )\n",
" _ = ax.loglog(\n",
" coupling_data[bench_name][in_index].x,\n",
" abs(coupling_data[bench_name][in_index].y),\n",
" label=\"Optickle\",\n",
" )\n",
" _ = ax.set_title(\"En phase\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" console.print()"
]
},
{
"cell_type": "markdown",
"id": "11606546-2606-404c-a75e-6e434d19084b",
"metadata": {},
"source": [
"### En fonction de la simulation"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "88002141-6a42-4ae0-a614-4c2b8dd336ac",
"metadata": {},
"outputs": [],
"source": [
"for bench_name in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
" in_index = 0\n",
" if bench_name == \"SDB1\":\n",
" in_index = 0\n",
" quad_index = (1 + in_index) % 2\n",
" Figure = figure(figsize=(14, 5))\n",
" _ = Figure.suptitle(\n",
" \"Comparaison des fonctions de transfert pour {}\".format(bench_name)\n",
" )\n",
" ax = Figure.add_subplot(1, 2, 1)\n",
" _ = ax.loglog(\n",
" quad_tf[bench_name].f,\n",
" abs(quad_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(quad_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"Quadrature de phase\",\n",
" )\n",
" _ = ax.loglog(\n",
" in_tf[bench_name].f,\n",
" abs(in_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(in_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"En phase\",\n",
" )\n",
" _ = ax.set_title(\"Finesse\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" ax = Figure.add_subplot(1, 2, 2)\n",
" _ = ax.loglog(\n",
" coupling_data[bench_name][quad_index].x,\n",
" abs(coupling_data[bench_name][quad_index].y),\n",
" label=\"Quadrature de phase\",\n",
" )\n",
" _ = ax.loglog(\n",
" coupling_data[bench_name][in_index].x,\n",
" abs(coupling_data[bench_name][in_index].y),\n",
" label=\"En phase\",\n",
" )\n",
" _ = ax.set_title(\"Optickle\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" show()"
]
},
{
"cell_type": "markdown",
"id": "0c70d12b-b5ae-44b7-b0d3-6f054b697300",
"metadata": {},
"source": [
"### En fonction du module/phase"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "cb1fd40c-83bb-4dc3-b259-7dd3a08d3b6b",
"metadata": {},
"outputs": [],
"source": [
"for bench_name in [\"SNEB\", \"SWEB\", \"SDB1\"]:\n",
" Figure = figure(figsize=(14, 5))\n",
" _ = Figure.suptitle(\n",
" \"Comparaison des fonctions de transfert pour {}\".format(bench_name)\n",
" )\n",
" ax = Figure.add_subplot(1, 2, 1)\n",
" _ = ax.loglog(\n",
" quad_tf[bench_name].f,\n",
" abs(quad_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(quad_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"Quadrature de phase\",\n",
" )\n",
" _ = ax.loglog(\n",
" in_tf[bench_name].f,\n",
" abs(in_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)])\n",
" / abs(in_tf[\"DARM\"][B1_detector, \"DARM\"])\n",
" / model.space_NI_NE.L.eval(),\n",
" label=\"En phase\",\n",
" )\n",
" _ = ax.set_title(\"Module\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" ax = Figure.add_subplot(1, 2, 2)\n",
" _ = ax.semilogx(\n",
" quad_tf[bench_name].f,\n",
" angle(quad_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)]) * 180 / pi,\n",
" label=\"Quadrature de phase\",\n",
" )\n",
" _ = ax.semilogx(\n",
" in_tf[bench_name].f,\n",
" angle(in_tf[bench_name][B1_detector, \"{}_z\".format(bench_name)]) * 180 / pi,\n",
" label=\"En phase\",\n",
" )\n",
" _ = ax.set_title(\"Finesse\")\n",
" _ = ax.legend()\n",
" ax.grid(True, \"both\", \"both\")\n",
" show()"
]
}
],
"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
}