""" Zernike Dashboard - Comprehensive Mirror Wavefront Analysis A complete, professional dashboard for M1 mirror optimization that combines: - Executive summary with all key metrics at a glance - All orientation views (40°, 60°, 90°) on a single page - MSF band analysis integrated - Light/white theme for better readability - Interactive comparison between subcases This is the unified replacement for the separate zernike_wfe and msf_zernike insights. Author: Atomizer Framework Date: 2024-12-22 """ from pathlib import Path from datetime import datetime from typing import Dict, Any, List, Optional, Tuple import numpy as np from math import factorial from numpy.linalg import LinAlgError from .base import StudyInsight, InsightConfig, InsightResult, register_insight # Lazy imports _plotly_loaded = False _go = None _make_subplots = None _Triangulation = None _OP2 = None _BDF = None _ZernikeOPDExtractor = None def _load_dependencies(): """Lazy load heavy dependencies.""" global _plotly_loaded, _go, _make_subplots, _Triangulation, _OP2, _BDF, _ZernikeOPDExtractor if not _plotly_loaded: import plotly.graph_objects as go from plotly.subplots import make_subplots from matplotlib.tri import Triangulation from pyNastran.op2.op2 import OP2 from pyNastran.bdf.bdf import BDF from ..extractors.extract_zernike_figure import ZernikeOPDExtractor _go = go _make_subplots = make_subplots _Triangulation = Triangulation _OP2 = OP2 _BDF = BDF _ZernikeOPDExtractor = ZernikeOPDExtractor _plotly_loaded = True # ============================================================================ # Color Themes # ============================================================================ LIGHT_THEME = { 'bg': '#ffffff', 'card_bg': '#f8fafc', 'card_border': '#e2e8f0', 'text': '#1e293b', 'text_muted': '#64748b', 'accent': '#3b82f6', 'success': '#10b981', 'warning': '#f59e0b', 'danger': '#ef4444', 'grid': 'rgba(100,116,139,0.2)', 'surface_bg': '#f1f5f9', } BAND_COLORS = { 'lsf': '#3b82f6', # Blue 'msf': '#f59e0b', # Amber 'hsf': '#ef4444', # Red } ANGLE_COLORS = { '40': '#10b981', # Emerald '60': '#3b82f6', # Blue '90': '#8b5cf6', # Purple } # ============================================================================ # Zernike Mathematics # ============================================================================ def noll_indices(j: int) -> Tuple[int, int]: """Convert Noll index to (n, m) radial/azimuthal orders.""" if j < 1: raise ValueError("Noll index j must be >= 1") count = 0 n = 0 while True: if n == 0: ms = [0] elif n % 2 == 0: ms = [0] + [m for k in range(1, n//2 + 1) for m in (-2*k, 2*k)] else: ms = [m for k in range(0, (n+1)//2) for m in (-(2*k+1), (2*k+1))] for m in ms: count += 1 if count == j: return n, m n += 1 def noll_to_radial_order(j: int) -> int: """Get radial order n for Noll index j.""" n, _ = noll_indices(j) return n def zernike_noll(j: int, r: np.ndarray, th: np.ndarray) -> np.ndarray: """Evaluate Zernike polynomial j at (r, theta).""" n, m = noll_indices(j) R = np.zeros_like(r) for s in range((n - abs(m)) // 2 + 1): c = ((-1)**s * factorial(n - s) / (factorial(s) * factorial((n + abs(m)) // 2 - s) * factorial((n - abs(m)) // 2 - s))) R += c * r**(n - 2*s) if m == 0: return R return R * (np.cos(m * th) if m > 0 else np.sin(-m * th)) def zernike_common_name(n: int, m: int) -> str: """Get common name for Zernike mode.""" names = { (0, 0): "Piston", (1, -1): "Tilt X", (1, 1): "Tilt Y", (2, 0): "Defocus", (2, -2): "Astig 45°", (2, 2): "Astig 0°", (3, -1): "Coma X", (3, 1): "Coma Y", (3, -3): "Trefoil X", (3, 3): "Trefoil Y", (4, 0): "Spherical", (4, -2): "Sec Astig X", (4, 2): "Sec Astig Y", (4, -4): "Quadrafoil X", (4, 4): "Quadrafoil Y", (5, -1): "Sec Coma X", (5, 1): "Sec Coma Y", (6, 0): "Sec Spherical", } return names.get((n, m), f"Z({n},{m})") def compute_zernike_coeffs( X: np.ndarray, Y: np.ndarray, vals: np.ndarray, n_modes: int, chunk_size: int = 100000 ) -> Tuple[np.ndarray, float]: """Fit Zernike coefficients to WFE data.""" Xc, Yc = X - np.mean(X), Y - np.mean(Y) R = float(np.max(np.hypot(Xc, Yc))) r = np.hypot(Xc / R, Yc / R).astype(np.float32) th = np.arctan2(Yc, Xc).astype(np.float32) mask = (r <= 1.0) & ~np.isnan(vals) if not np.any(mask): raise RuntimeError("No valid points inside unit disk.") idx = np.nonzero(mask)[0] m = int(n_modes) G = np.zeros((m, m), dtype=np.float64) h = np.zeros((m,), dtype=np.float64) v = vals.astype(np.float64) for start in range(0, len(idx), chunk_size): sl = idx[start:start + chunk_size] r_b, th_b, v_b = r[sl], th[sl], v[sl] Zb = np.column_stack([zernike_noll(j, r_b, th_b).astype(np.float32) for j in range(1, m + 1)]) G += (Zb.T @ Zb).astype(np.float64) h += (Zb.T @ v_b).astype(np.float64) try: coeffs = np.linalg.solve(G, h) except LinAlgError: coeffs = np.linalg.lstsq(G, h, rcond=None)[0] return coeffs, R def classify_modes_by_band( n_modes: int, lsf_max: int = 10, msf_max: int = 50 ) -> Dict[str, List[int]]: """Classify Zernike modes into LSF/MSF/HSF bands.""" bands = {'lsf': [], 'msf': [], 'hsf': []} for j in range(1, n_modes + 1): n = noll_to_radial_order(j) if n <= lsf_max: bands['lsf'].append(j) elif n <= msf_max: bands['msf'].append(j) else: bands['hsf'].append(j) return bands def compute_band_rss(coeffs: np.ndarray, bands: Dict[str, List[int]]) -> Dict[str, float]: """Compute RSS (Root Sum Square) of coefficients in each band.""" rss = {} for band_name, indices in bands.items(): band_coeffs = [coeffs[j-1] for j in indices if j-1 < len(coeffs)] rss[band_name] = float(np.sqrt(np.sum(np.array(band_coeffs)**2))) rss['total'] = float(np.sqrt(np.sum(coeffs**2))) # Filtered: exclude J1-J4 (piston, tip, tilt, defocus) filtered_coeffs = coeffs[4:] if len(coeffs) > 4 else np.array([]) rss['filtered'] = float(np.sqrt(np.sum(filtered_coeffs**2))) return rss # ============================================================================ # Default Configuration # ============================================================================ DEFAULT_CONFIG = { 'n_modes': 50, 'amp': 0.5, 'pancake': 3.0, 'plot_downsample': 8000, 'filter_low_orders': 4, 'lsf_max': 10, 'msf_max': 50, 'disp_unit': 'mm', 'use_opd': True, # Use OPD method by default } @register_insight class ZernikeDashboardInsight(StudyInsight): """ Comprehensive Zernike Dashboard for M1 Mirror Optimization. Generates a single-page dashboard with: - Executive summary with key metrics - All orientations (40°, 60°, 90°) in comparison view - MSF band analysis - Interactive surface plots - Professional light theme """ insight_type = "zernike_dashboard" name = "Zernike Dashboard" description = "Comprehensive mirror WFE dashboard with all orientations and MSF analysis" category = "optical" applicable_to = ["mirror", "optics", "wfe"] required_files = ["*.op2"] def __init__(self, study_path: Path): super().__init__(study_path) self.op2_path: Optional[Path] = None self.geo_path: Optional[Path] = None self._node_geo: Optional[Dict] = None self._displacements: Optional[Dict] = None def can_generate(self) -> bool: """Check if OP2 and geometry files exist.""" search_paths = [ self.results_path, self.study_path / "2_iterations", self.setup_path / "model", ] for search_path in search_paths: if not search_path.exists(): continue op2_files = list(search_path.glob("**/*solution*.op2")) if not op2_files: op2_files = list(search_path.glob("**/*.op2")) if op2_files: self.op2_path = max(op2_files, key=lambda p: p.stat().st_mtime) break if self.op2_path is None: return False try: self.geo_path = self._find_geometry_file(self.op2_path) return True except FileNotFoundError: return False def _find_geometry_file(self, op2_path: Path) -> Path: """Find BDF/DAT geometry file for OP2.""" folder = op2_path.parent base = op2_path.stem for ext in ['.dat', '.bdf']: cand = folder / (base + ext) if cand.exists(): return cand for f in folder.iterdir(): if f.suffix.lower() in ['.dat', '.bdf']: return f raise FileNotFoundError(f"No geometry file found for {op2_path}") def _load_data(self): """Load geometry and displacement data.""" if self._node_geo is not None: return _load_dependencies() bdf = _BDF() bdf.read_bdf(str(self.geo_path)) self._node_geo = {int(nid): node.get_position() for nid, node in bdf.nodes.items()} op2 = _OP2() op2.read_op2(str(self.op2_path)) if not op2.displacements: raise RuntimeError("No displacement data in OP2") self._displacements = {} for key, darr in op2.displacements.items(): data = darr.data dmat = data[0] if data.ndim == 3 else (data if data.ndim == 2 else None) if dmat is None: continue ngt = darr.node_gridtype.astype(int) node_ids = ngt if ngt.ndim == 1 else ngt[:, 0] isubcase = getattr(darr, 'isubcase', None) label = str(isubcase) if isubcase else str(key) self._displacements[label] = { 'node_ids': node_ids.astype(int), 'disp': dmat.copy() } def _build_wfe_arrays_standard( self, label: str, disp_unit: str = 'mm' ) -> Dict[str, np.ndarray]: """Build arrays using standard Z-only method.""" nm_per_unit = 1e6 if disp_unit == 'mm' else 1e9 data = self._displacements[label] node_ids = data['node_ids'] dmat = data['disp'] X, Y, WFE = [], [], [] valid_nids = [] dx_arr, dy_arr, dz_arr = [], [], [] for nid, vec in zip(node_ids, dmat): geo = self._node_geo.get(int(nid)) if geo is None: continue X.append(geo[0]) Y.append(geo[1]) wfe = vec[2] * 2.0 * nm_per_unit WFE.append(wfe) valid_nids.append(nid) dx_arr.append(vec[0]) dy_arr.append(vec[1]) dz_arr.append(vec[2]) return { 'X': np.array(X), 'Y': np.array(Y), 'WFE': np.array(WFE), 'node_ids': np.array(valid_nids), 'dx': np.array(dx_arr), 'dy': np.array(dy_arr), 'dz': np.array(dz_arr), 'lateral_disp': np.sqrt(np.array(dx_arr)**2 + np.array(dy_arr)**2), } def _build_wfe_arrays_opd( self, label: str, disp_unit: str = 'mm' ) -> Dict[str, np.ndarray]: """Build arrays using OPD method (recommended).""" _load_dependencies() extractor = _ZernikeOPDExtractor( self.op2_path, bdf_path=self.geo_path, displacement_unit=disp_unit ) opd_data = extractor._build_figure_opd_data(label) return { 'X': opd_data['x_deformed'], 'Y': opd_data['y_deformed'], 'WFE': opd_data['wfe_nm'], 'node_ids': opd_data['node_ids'], 'dx': opd_data['dx'], 'dy': opd_data['dy'], 'dz': opd_data['dz'], 'lateral_disp': opd_data['lateral_disp'], 'x_original': opd_data['x_original'], 'y_original': opd_data['y_original'], } def _build_wfe_arrays( self, label: str, disp_unit: str = 'mm', use_opd: bool = True ) -> Dict[str, np.ndarray]: """Build X, Y, WFE arrays for a subcase.""" if use_opd: return self._build_wfe_arrays_opd(label, disp_unit) else: return self._build_wfe_arrays_standard(label, disp_unit) def _compute_relative_wfe( self, data1: Dict[str, np.ndarray], data2: Dict[str, np.ndarray] ) -> Dict[str, np.ndarray]: """Compute relative WFE (target - reference) for common nodes.""" X1, Y1, WFE1, nids1 = data1['X'], data1['Y'], data1['WFE'], data1['node_ids'] X2, Y2, WFE2, nids2 = data2['X'], data2['Y'], data2['WFE'], data2['node_ids'] ref_map = {int(nid): (x, y, w) for nid, x, y, w in zip(nids2, X2, Y2, WFE2)} dx1_map = {int(nid): d for nid, d in zip(data1['node_ids'], data1['dx'])} dy1_map = {int(nid): d for nid, d in zip(data1['node_ids'], data1['dy'])} dz1_map = {int(nid): d for nid, d in zip(data1['node_ids'], data1['dz'])} dx2_map = {int(nid): d for nid, d in zip(data2['node_ids'], data2['dx'])} dy2_map = {int(nid): d for nid, d in zip(data2['node_ids'], data2['dy'])} dz2_map = {int(nid): d for nid, d in zip(data2['node_ids'], data2['dz'])} X_rel, Y_rel, WFE_rel, nids_rel = [], [], [], [] dx_rel, dy_rel, dz_rel = [], [], [] for nid, x, y, w in zip(nids1, X1, Y1, WFE1): nid = int(nid) if nid not in ref_map or nid not in dx2_map: continue _, _, w_ref = ref_map[nid] X_rel.append(x) Y_rel.append(y) WFE_rel.append(w - w_ref) nids_rel.append(nid) dx_rel.append(dx1_map[nid] - dx2_map[nid]) dy_rel.append(dy1_map[nid] - dy2_map[nid]) dz_rel.append(dz1_map[nid] - dz2_map[nid]) return { 'X': np.array(X_rel), 'Y': np.array(Y_rel), 'WFE': np.array(WFE_rel), 'node_ids': np.array(nids_rel), 'dx': np.array(dx_rel), 'dy': np.array(dy_rel), 'dz': np.array(dz_rel), 'lateral_disp': np.sqrt(np.array(dx_rel)**2 + np.array(dy_rel)**2) if dx_rel else np.array([]), } def _compute_metrics( self, X: np.ndarray, Y: np.ndarray, W_nm: np.ndarray, n_modes: int, filter_orders: int, cfg: Dict ) -> Dict[str, Any]: """Compute comprehensive metrics including band analysis.""" coeffs, R = compute_zernike_coeffs(X, Y, W_nm, n_modes) Xc = X - np.mean(X) Yc = Y - np.mean(Y) r = np.hypot(Xc / R, Yc / R) th = np.arctan2(Yc, Xc) Z = np.column_stack([zernike_noll(j, r, th) for j in range(1, n_modes + 1)]) W_res_filt = W_nm - Z[:, :filter_orders].dot(coeffs[:filter_orders]) W_res_filt_j1to3 = W_nm - Z[:, :3].dot(coeffs[:3]) # Band analysis bands = classify_modes_by_band(n_modes, cfg['lsf_max'], cfg['msf_max']) band_rss = compute_band_rss(coeffs, bands) # Aberration magnitudes aberrations = { 'defocus': float(abs(coeffs[3])) if len(coeffs) > 3 else 0.0, 'astigmatism': float(np.sqrt(coeffs[4]**2 + coeffs[5]**2)) if len(coeffs) > 5 else 0.0, 'coma': float(np.sqrt(coeffs[6]**2 + coeffs[7]**2)) if len(coeffs) > 7 else 0.0, 'trefoil': float(np.sqrt(coeffs[8]**2 + coeffs[9]**2)) if len(coeffs) > 9 else 0.0, 'spherical': float(abs(coeffs[10])) if len(coeffs) > 10 else 0.0, } return { 'coefficients': coeffs, 'R': R, 'global_rms': float(np.sqrt(np.mean(W_nm**2))), 'filtered_rms': float(np.sqrt(np.mean(W_res_filt**2))), 'rms_j1to3': float(np.sqrt(np.mean(W_res_filt_j1to3**2))), 'W_res_filt': W_res_filt, 'band_rss': band_rss, 'aberrations': aberrations, } def _generate_html( self, all_data: Dict[str, Dict], cfg: Dict, timestamp: str ) -> str: """Generate comprehensive HTML dashboard.""" _load_dependencies() theme = LIGHT_THEME # Extract key metrics for executive summary metrics_40 = all_data.get('40_vs_20', {}).get('metrics', {}) metrics_60 = all_data.get('60_vs_20', {}).get('metrics', {}) metrics_90 = all_data.get('90_abs', {}).get('metrics', {}) rms_40 = metrics_40.get('filtered_rms', 0) rms_60 = metrics_60.get('filtered_rms', 0) rms_90_mfg = metrics_90.get('rms_j1to3', 0) # Generate HTML html = f''' Zernike Dashboard - Atomizer

Zernike WFE Dashboard

M1 Mirror Optimization Analysis • OPD Method (X,Y,Z corrected)
Generated: {timestamp}
40° vs 20° RMS
{rms_40:.2f} nm
Target: ≤ 4.0 nm
60° vs 20° RMS
{rms_60:.2f} nm
Target: ≤ 10.0 nm
90° MFG (J1-J3)
{rms_90_mfg:.2f} nm
Target: ≤ 20.0 nm
Weighted Score
{self._compute_weighted_score(rms_40, rms_60, rms_90_mfg):.1f}
Lower is better
Detailed Metrics Comparison
All values in nm • J1-J4 filtered (piston, tip, tilt, defocus removed)
Metric 40° vs 20° 60° vs 20° 90° MFG Target
Filtered RMS (J1-J4) {rms_40:.2f} {rms_60:.2f} ≤ 4 / 10 nm
Optician Workload (J1-J3) {rms_90_mfg:.2f} ≤ 20 nm
Astigmatism (J5+J6) {metrics_40.get('aberrations', {}).get('astigmatism', 0):.2f} {metrics_60.get('aberrations', {}).get('astigmatism', 0):.2f} {metrics_90.get('aberrations', {}).get('astigmatism', 0):.2f}
Coma (J7+J8) {metrics_40.get('aberrations', {}).get('coma', 0):.2f} {metrics_60.get('aberrations', {}).get('coma', 0):.2f} {metrics_90.get('aberrations', {}).get('coma', 0):.2f}
Trefoil (J9+J10) {metrics_40.get('aberrations', {}).get('trefoil', 0):.2f} {metrics_60.get('aberrations', {}).get('trefoil', 0):.2f} {metrics_90.get('aberrations', {}).get('trefoil', 0):.2f}
Spherical (J11) {metrics_40.get('aberrations', {}).get('spherical', 0):.2f} {metrics_60.get('aberrations', {}).get('spherical', 0):.2f} {metrics_90.get('aberrations', {}).get('spherical', 0):.2f}
Spatial Frequency Band Analysis
RSS decomposition • LSF: n≤10 (M2 correctable) • MSF: n=11-50 (support print-through) • HSF: n>50
Band 40° vs 20° 60° vs 20° 90° MFG Physics
LSF (n ≤ 10) {metrics_40.get('band_rss', {}).get('lsf', 0):.2f} nm {metrics_60.get('band_rss', {}).get('lsf', 0):.2f} nm {metrics_90.get('band_rss', {}).get('lsf', 0):.2f} nm M2 hexapod correctable
MSF (n = 11-50) {metrics_40.get('band_rss', {}).get('msf', 0):.2f} nm {metrics_60.get('band_rss', {}).get('msf', 0):.2f} nm {metrics_90.get('band_rss', {}).get('msf', 0):.2f} nm Support print-through
HSF (n > 50) {metrics_40.get('band_rss', {}).get('hsf', 0):.2f} nm {metrics_60.get('band_rss', {}).get('hsf', 0):.2f} nm {metrics_90.get('band_rss', {}).get('hsf', 0):.2f} nm High frequency (mesh limit)
Total RSS {metrics_40.get('band_rss', {}).get('total', 0):.2f} nm {metrics_60.get('band_rss', {}).get('total', 0):.2f} nm {metrics_90.get('band_rss', {}).get('total', 0):.2f} nm
Surface Residual Visualization
WFE residual after J1-J4 removal • Visual amplification applied for clarity
40° vs 20° (Relative)
60° vs 20° (Relative)
90° Manufacturing (Absolute)
Zernike Coefficient Magnitudes
Absolute values • Color-coded by spatial frequency band
LSF (n ≤ 10)
MSF (n = 11-50)
HSF (n > 50)
''' return html def _compute_weighted_score(self, rms_40: float, rms_60: float, rms_90: float) -> float: """Compute weighted optimization score.""" # Standard weights from M1 optimization w40 = 100 / 4.0 # Target 4 nm w60 = 50 / 10.0 # Target 10 nm w90 = 20 / 20.0 # Target 20 nm return w40 * rms_40 + w60 * rms_60 + w90 * rms_90 def _get_status_class(self, value: float, target: float, warn: float) -> str: """Get CSS class for value status.""" if value <= target: return ' pass' elif value <= warn: return ' warn' else: return ' fail' def _generate_plot_js(self, all_data: Dict[str, Dict], cfg: Dict) -> str: """Generate JavaScript for Plotly plots.""" _load_dependencies() amp = cfg.get('amp', 0.5) downsample = cfg.get('plot_downsample', 8000) js_parts = [] # Generate 3D surface plots for each angle for key, angle, div_id in [ ('40_vs_20', '40', 'plot-40'), ('60_vs_20', '60', 'plot-60'), ('90_abs', '90', 'plot-90'), ]: data = all_data.get(key, {}) X = data.get('X', np.array([])) Y = data.get('Y', np.array([])) W = data.get('metrics', {}).get('W_res_filt', np.array([])) if len(X) == 0: continue # Downsample n = len(X) if n > downsample: rng = np.random.default_rng(42) sel = rng.choice(n, size=downsample, replace=False) Xp, Yp, Wp = X[sel], Y[sel], W[sel] else: Xp, Yp, Wp = X, Y, W res_amp = amp * Wp max_amp = float(np.max(np.abs(res_amp))) if res_amp.size else 1.0 # Build triangulation try: tri = _Triangulation(Xp, Yp) i_arr = tri.triangles[:, 0].tolist() j_arr = tri.triangles[:, 1].tolist() k_arr = tri.triangles[:, 2].tolist() except: i_arr, j_arr, k_arr = [], [], [] color = ANGLE_COLORS[angle] js_parts.append(f''' // {angle}° Plot Plotly.newPlot('{div_id}', [{{ type: 'mesh3d', x: {Xp.tolist()}, y: {Yp.tolist()}, z: {res_amp.tolist()}, i: {i_arr}, j: {j_arr}, k: {k_arr}, intensity: {res_amp.tolist()}, colorscale: 'RdBu_r', showscale: true, colorbar: {{ title: 'WFE (nm)', thickness: 12, len: 0.6 }}, lighting: {{ ambient: 0.5, diffuse: 0.8, specular: 0.3 }}, hovertemplate: 'X: %{{x:.1f}} mm
Y: %{{y:.1f}} mm
WFE: %{{z:.2f}} nm' }}], {{ margin: {{t: 10, b: 10, l: 10, r: 10}}, scene: {{ camera: {{eye: {{x: 1.2, y: 1.2, z: 0.8}}}}, xaxis: {{title: 'X (mm)', showgrid: true, gridcolor: 'rgba(100,100,100,0.2)'}}, yaxis: {{title: 'Y (mm)', showgrid: true, gridcolor: 'rgba(100,100,100,0.2)'}}, zaxis: {{title: 'WFE (nm)', range: [{-max_amp * 3}, {max_amp * 3}], showgrid: true, gridcolor: 'rgba(100,100,100,0.2)'}}, aspectmode: 'manual', aspectratio: {{x: 1, y: 1, z: 0.4}}, bgcolor: '#f8fafc' }}, paper_bgcolor: '#ffffff' }}, {{responsive: true}}); ''') # Generate coefficient bar chart metrics_40 = all_data.get('40_vs_20', {}).get('metrics', {}) coeffs = metrics_40.get('coefficients', np.array([])) n_modes = len(coeffs) if len(coeffs) > 0 else 50 if len(coeffs) > 0: bands = classify_modes_by_band(n_modes, cfg.get('lsf_max', 10), cfg.get('msf_max', 50)) bar_colors = [] labels = [] for j in range(1, n_modes + 1): n, m = noll_indices(j) labels.append(f'J{j} {zernike_common_name(n, m)}') if j in bands['lsf']: bar_colors.append(BAND_COLORS['lsf']) elif j in bands['msf']: bar_colors.append(BAND_COLORS['msf']) else: bar_colors.append(BAND_COLORS['hsf']) coeff_abs = np.abs(coeffs).tolist() js_parts.append(f''' // Coefficient Bar Chart Plotly.newPlot('coeff-chart', [{{ type: 'bar', x: {coeff_abs}, y: {labels}, orientation: 'h', marker: {{ color: {bar_colors} }}, hovertemplate: '%{{y}}
|Coeff| = %{{x:.3f}} nm' }}], {{ margin: {{t: 20, b: 40, l: 200, r: 20}}, height: 800, xaxis: {{title: '|Coefficient| (nm)', gridcolor: 'rgba(100,100,100,0.2)'}}, yaxis: {{autorange: 'reversed'}}, paper_bgcolor: '#ffffff', plot_bgcolor: '#f8fafc' }}, {{responsive: true}}); ''') return '\n'.join(js_parts) def _generate(self, config: InsightConfig) -> InsightResult: """Generate comprehensive Zernike dashboard.""" self._load_data() cfg = {**DEFAULT_CONFIG, **config.extra} n_modes = cfg['n_modes'] filter_orders = cfg['filter_low_orders'] disp_unit = cfg['disp_unit'] use_opd = cfg['use_opd'] # Map subcases disps = self._displacements if '1' in disps and '2' in disps: sc_map = {'90': '1', '20': '2', '40': '3', '60': '4'} elif '90' in disps and '20' in disps: sc_map = {'90': '90', '20': '20', '40': '40', '60': '60'} else: available = sorted(disps.keys(), key=lambda x: int(x) if x.isdigit() else 0) if len(available) >= 4: sc_map = {'90': available[0], '20': available[1], '40': available[2], '60': available[3]} else: return InsightResult(success=False, error=f"Need 4 subcases, found: {available}") # Validate for angle, label in sc_map.items(): if label not in disps: return InsightResult(success=False, error=f"Subcase '{label}' (angle {angle}) not found") timestamp = datetime.now().strftime("%Y-%m-%d %H:%M:%S") output_dir = config.output_dir or self.insights_path output_dir.mkdir(parents=True, exist_ok=True) all_data = {} # Load reference data (20°) data_ref = self._build_wfe_arrays(sc_map['20'], disp_unit, use_opd) # Process each angle for angle in ['40', '60']: data = self._build_wfe_arrays(sc_map[angle], disp_unit, use_opd) rel_data = self._compute_relative_wfe(data, data_ref) metrics = self._compute_metrics( rel_data['X'], rel_data['Y'], rel_data['WFE'], n_modes, filter_orders, cfg ) all_data[f'{angle}_vs_20'] = { 'X': rel_data['X'], 'Y': rel_data['Y'], 'WFE': rel_data['WFE'], 'metrics': metrics, } # 90° absolute (manufacturing) data_90 = self._build_wfe_arrays(sc_map['90'], disp_unit, use_opd) metrics_90 = self._compute_metrics( data_90['X'], data_90['Y'], data_90['WFE'], n_modes, filter_orders, cfg ) all_data['90_abs'] = { 'X': data_90['X'], 'Y': data_90['Y'], 'WFE': data_90['WFE'], 'metrics': metrics_90, } # Generate HTML html = self._generate_html(all_data, cfg, timestamp) file_timestamp = datetime.now().strftime("%Y%m%d_%H%M%S") html_path = output_dir / f"zernike_dashboard_{file_timestamp}.html" html_path.write_text(html, encoding='utf-8') # Build summary summary = { '40_vs_20_filtered_rms': all_data['40_vs_20']['metrics']['filtered_rms'], '60_vs_20_filtered_rms': all_data['60_vs_20']['metrics']['filtered_rms'], '90_mfg_rms_j1to3': all_data['90_abs']['metrics']['rms_j1to3'], '90_mfg_filtered_rms': all_data['90_abs']['metrics']['filtered_rms'], 'weighted_score': self._compute_weighted_score( all_data['40_vs_20']['metrics']['filtered_rms'], all_data['60_vs_20']['metrics']['filtered_rms'], all_data['90_abs']['metrics']['rms_j1to3'] ), 'method': 'OPD' if use_opd else 'Standard', 'n_modes': n_modes, } return InsightResult( success=True, html_path=html_path, summary=summary )