function [fitResults, rawCurves] = fitTwoGaussianCurves(S_theta_all, theta_vals, varargin) %% fitTwoGaussianCurves % Fits a two-Gaussian model to multiple spectral curves from θ=0 up to % the first secondary peak ≥50% of the primary peak, within 0–π/2 radians. % % Author: Karthik % Date: 2025-10-06 % Version: 1.1 % % Inputs: % S_theta_all - 1xN cell array of spectral curves % theta_vals - vector of corresponding theta values (radians) % % Optional name-value pairs: % 'MaxTheta' - Maximum theta for search (default: pi/2) % 'ResidualThreshold' - Max relative deviation allowed (default: 0.3) % 'PositionThreshold' - Minimum angular separation between main and % secondary peaks (default: 0.1 rad) % 'AmplitudeThreshold' - Minimum fraction of primary peak amplitude for secondary peak (default: 0.5) % % Output: % fitResults - struct array with fields: % .pFit, .thetaFit, .xFit, .yFit, .thetaFine, .isValid, .fitMaxTheta % --- Parse optional inputs --- p = inputParser; addParameter(p, 'MaxTheta', pi/2, @(x) isnumeric(x) && isscalar(x)); addParameter(p, 'ResidualThreshold', 0.3, @(x) isnumeric(x) && isscalar(x)); addParameter(p, 'PositionThreshold', 0.1, @(x) isnumeric(x) && isscalar(x)); addParameter(p, 'AmplitudeThreshold', 0.5, @(x) isnumeric(x) && isscalar(x)); addParameter(p, 'RecenterCurves', true, @islogical); parse(p, varargin{:}); opts = p.Results; Ncurves = numel(S_theta_all); fitResults = struct('pFit',[],'thetaFit',[],'xFit',[],'yFit',[],... 'thetaFine',[],'isValid',[],'fitMaxTheta',[]); rawCurves = struct('x', cell(1, Ncurves), 'theta', cell(1, Ncurves)); if opts.RecenterCurves % --- Preprocess curves: shift first peak to zero without wrapping --- S_theta_all_shifted = cell(size(S_theta_all)); for k = 1:Ncurves curve = S_theta_all{k}; idx_range = find(theta_vals >= 0 & theta_vals <= opts.MaxTheta); if isempty(idx_range) [~, peak_idx] = max(curve); else [~, local_max_idx] = max(curve(idx_range)); peak_idx = idx_range(local_max_idx); end S_theta_all_shifted{k} = curve(peak_idx:end); end % --- Pad curves to same length --- Npoints_shifted = max(cellfun(@numel, S_theta_all_shifted)); for k = 1:Ncurves len = numel(S_theta_all_shifted{k}); if len < Npoints_shifted S_theta_all_shifted{k} = [S_theta_all_shifted{k}, nan(1, Npoints_shifted-len)]; end end theta_recentered = theta_vals(1:Npoints_shifted) - theta_vals(1); else S_theta_all_shifted = S_theta_all; theta_recentered = theta_vals; end for k = 1:Ncurves x = S_theta_all_shifted{k}; theta = theta_recentered; validIdx = ~isnan(x); x = x(validIdx); theta = theta(validIdx); mask = theta>=0 & theta<=opts.MaxTheta; x = x(mask); theta = theta(mask); x = x / x(1); % --- Store raw data for plotting --- rawCurves(k).x = x; rawCurves(k).theta = theta; try % --- Filter 1: Ensure curve covers full 0–MaxTheta range --- thetaRange = max(theta) - min(theta); if thetaRange < 0.9 * opts.MaxTheta warning('Curve %d excluded: does not span full 0–MaxTheta range.', k); fitResults(k) = fillZeroStruct(theta, opts.MaxTheta); continue; end if numel(theta) < 8 warning('Curve %d too short (<8 points), skipping.', k); fitResults(k) = fillNaNStruct(); continue; end xSmooth = smooth(x, 5, 'moving'); curveAmp = max(xSmooth) - min(xSmooth); %% 1) Attempt primary detection via findpeaks minProm = opts.AmplitudeThreshold * curveAmp; [pk, locIdx] = findpeaks(xSmooth, 'MinPeakProminence', minProm); thetaPeaks = theta(locIdx); foundSecondary = false; fitMaxTheta = []; secondaryAmp = []; secondaryIdx = find(thetaPeaks > opts.PositionThreshold, 1, 'first'); if ~isempty(secondaryIdx) foundSecondary = true; fitMaxTheta = thetaPeaks(secondaryIdx); secondaryAmp = pk(secondaryIdx); end %% 2) If secondary found, fit constrained two-Gaussian directly if foundSecondary twoGauss = @(p,theta) ... 1.0 * exp(-0.5*((theta - p(1))/max(p(2),1e-6)).^2) + ... p(3) * exp(-0.5*((theta - p(4))/max(p(5),1e-6)).^2); mu1_guess = 0; sigma1_guess = 0.1 * fitMaxTheta; A2_guess = secondaryAmp; mu2_guess = fitMaxTheta; sigma2_guess = fitMaxTheta / 3; p0 = [mu1_guess, sigma1_guess, A2_guess, mu2_guess, sigma2_guess]; lb = [0, 1e-6, 0, max(fitMaxTheta - 0.1, 0), 1e-6]; ub = [fitMaxTheta/4, fitMaxTheta, 2, min(fitMaxTheta + 0.1, opts.MaxTheta), fitMaxTheta]; optsLSQ = optimoptions('lsqcurvefit','Display','off', ... 'MaxFunctionEvaluations',1e4,'MaxIterations',1e4); try pFit = lsqcurvefit(twoGauss, p0, theta, x, lb, ub, optsLSQ); % --- Handle missing or weak secondary peak, discard poor fits --- if pFit(3) <= opts.AmplitudeThreshold % amplitude too small → discard secondary pFit(3) = 0; pFit(4:5) = NaN; else % --- Compute RMS of fit vs raw curve in [0, pi/4] --- maskRMS = theta >= 0 & theta <= pi/4; if any(maskRMS) x_raw_subset = x(maskRMS); y_fit_subset = twoGauss(pFit, theta(maskRMS)); rms_val = sqrt(mean((x_raw_subset - y_fit_subset).^2)); if rms_val > opts.ResidualThreshold % residual too large → discard secondary pFit(3) = NaN; pFit(4:5) = NaN; end end end % --- Map to canonical format with primary amplitude = 1 --- pCanonical = [1.0, pFit(1), pFit(2), pFit(3), pFit(4), pFit(5)]; thetaFine = linspace(0, opts.MaxTheta, 500); yFit = 1.0*exp(-0.5*((thetaFine - pCanonical(2))/max(pCanonical(3),1e-6)).^2) + ... pCanonical(4)*exp(-0.5*((thetaFine - pCanonical(5))/max(pCanonical(6),1e-6)).^2); fitResults(k).pFit = pCanonical; fitResults(k).thetaFit = theta; fitResults(k).xFit = x; fitResults(k).thetaFine = thetaFine; fitResults(k).yFit = yFit; fitResults(k).isValid = true; fitResults(k).fitMaxTheta = fitMaxTheta; continue; % success → next curve catch warning('Curve %d constrained fit failed, switching to generic fit.', k); end end %% 3) Generic two-Gaussian fallback (no or failed secondary detection) try twoGaussFree = @(p,th) ... p(1)*exp(-0.5*((th - p(2))/max(p(3),1e-6)).^2) + ... p(4)*exp(-0.5*((th - p(5))/max(p(6),1e-6)).^2); p0 = [1, 0, 0.1, 0.3, opts.MaxTheta/3, 0.1]; lb = [0, 0, 1e-6, 0, opts.PositionThreshold, 1e-6]; ub = [2, opts.MaxTheta/2, opts.MaxTheta, 2, opts.MaxTheta, opts.MaxTheta]; optsLSQfree = optimoptions('lsqcurvefit', 'Display', 'off', ... 'MaxFunctionEvaluations', 2e4, 'MaxIterations', 2e4); pFree = lsqcurvefit(twoGaussFree, p0, theta(:), xSmooth(:), lb, ub, optsLSQfree); % --- Handle missing secondary peak --- if pFree(4) <= 1e-3 pFree(4) = 0; % secondary amplitude = 0 pFree(5:6) = NaN; % secondary position & width = NaN end fitMaxTheta = max(pFree(5), pFree(2)); thetaFine = linspace(0, opts.MaxTheta, 500); yFit = twoGaussFree(pFree, thetaFine); fitResults(k).pFit = pFree; fitResults(k).thetaFit = theta; fitResults(k).xFit = x; fitResults(k).thetaFine = thetaFine; fitResults(k).yFit = yFit; fitResults(k).isValid = true; fitResults(k).fitMaxTheta = fitMaxTheta; catch warning('Curve %d: generic two-Gaussian fit failed.', k); fitResults(k) = fillNaNStruct(); end catch warning('Curve %d fit failed completely.', k); fitResults(k) = fillNaNStruct(); end end end % --- Helper: return struct with NaNs (fit failed) --- function s = fillNaNStruct() s = struct( ... 'pFit', nan(1,5), ... 'thetaFit', nan, ... 'xFit', nan, ... 'yFit', nan, ... 'thetaFine', nan, ... 'isValid', false, ... 'fitMaxTheta', nan ... ); end % --- Helper: return struct with zeros (no secondary peak found) --- function s = fillZeroStruct(theta, fitMaxTheta) s = struct( ... 'pFit', zeros(1,5), ... 'thetaFit', theta, ... 'xFit', zeros(size(theta)), ... 'yFit', zeros(size(theta)), ... 'thetaFine', zeros(1,500), ... 'isValid', false, ... 'fitMaxTheta', fitMaxTheta ... ); end