function varargout=swdisk(m,N,K,L,x,method,scalem) % [E,V,Nm,c,x,E2]=SWDISK(m,N,K,L,x,method,scalem) % % Calculates the radial part of the Slepian functions concentrated to a % circular domain in Cartesian space. For comparison, can solve the % integral equation directly, and we can also compare with fully % two-dimensional methods using LOCALIZATION2D. Note that increasing L is % almost always a good idea. Remember the distinction between the overall % Shannon number N and the fixed-order Shannon number Nm, and see NSUBM. % % INPUT: % % m The fixed order of the eigenvalue problem % N The Shannon number of the full concentration problem [default: 20] % K The number of functions (best first) that you want [default: N] % L The maximum degree in the Jacobi/Bessel expansion [default: 2N] % or else the maximum polynomial degree for which GL will be exact % x The abscissas required [default: 0-1 for 'DV' and 0-2 for 'SE'] % If NaN only computes the eigenvalues % method 'DV' by the method of De Villiers % 'SE' by Slepian "extension" [default] % 'GL' by direct Gauss-Legendre integration % scalem 1 Scales the solution for weightless orthogonality [default] % 0 Scales the solution for orthogonality with weight x % % OUTPUT: % % E The orthogonal radial Slepian eigenfunctions for the disk, % where 2\pi\int_0^1 E_i(x) E_j(x) \,dx=\delta_{ij}V_i OR % 2\pi\int_0^1 E_i(x) E_j(x) x\,dx=\delta_{ij}V_in % V The concentration eigenvalues % N The partial Shannon number for this order % c The concentration factor c=2*sqrt(N) % x The abscissas at which the functions are evaluated % E2 The orthonormal radial Slepian eigenfunctions for the disk, % where \int_0^1 E_i(x) E_j(x) \,dx=\delta_{ij} OR % \int_0^1 E_i(x) E_j(x) x\,dx=\delta_{ij} % % SEE ALSO: % % JACOBI, LOCALIZATION2D, KERNELC2D, DEVILLIERS % % EXAMPLE: % % swdisk('demoX') where X=1->10 % % Last modified by dongwang-at-princeton.edu, 08/06/2008 % Last modified by fjsimons-at-alum.mit.edu, 03/03/2011 % Supply defaults defval('m',3) defval('N',3) defval('K',ceil(N)) defval('method','SE') defval('scalem',1) % Some heuristic defaults - make identical in DEVILLIERS defval('L',min(max([2*K ceil(2*N) 2*m 10]),... 84+192*sum([method=='GL']))) % Time the code t0=clock; if ~isstr(m) switch method case 'DV' disp('Jacobi expansion method on the interval') % Get the Jacobi expansion coefficients [D,V,g,c]=devilliers(m,N,K,L); % Quote the maximum degree of the Jacobi expansion disp(sprintf('\nJacobi expansion up to degree %i\n',L)) % What is the interval over which we are looking? defval('x',linspace(0,1,200)); if ~isnan(x) if any(x<0) error('Arguments must be zero or positive') end if ~~sum(x(:)>1 | x(:)<0); warning('Results will be inaccurate outside [0 1]'); end % Calculate the required set of scaled Jacobi polynomials... % ... for the normalization... [w,xGL,NGL]=gausslegendrecof(L,[],[0 1]); PGL=getscaledjacobis(xGL,m,L); % Calculate the required set of scaled Jacobi polynomials... % ... for the actual points requested P=getscaledjacobis(x,m,L); % Calculate the whole set functions as a matrix with row dimension the % evaluation points and column dimension the eigenfunction rank. EGL=PGL'*D; E=P'*D; % Calculate the normalization on the unit interval orv2=2*pi*EGL'*diag(w)*EGL; disp(sprintf('\nGauss-Legendre integration, Nystrom method with %i nodes',... NGL)) disp(sprintf('Normalization: Mean absolute error %8.3e',... mean(mean(abs(orv2-diag(diag(orv2))))))) % Sometimes you won't want this, and undoing might be hard when the % eigenvalues are very small, they may be inaccurate E2=sqrt(2*pi)*[diag(1./sqrt(diag(orv2)))*E']'; % Apply and return the normalized result E=[diag(sqrt(V(:))./sqrt(diag(orv2)))*E']'; else E=NaN; E2=NaN; end case 'SE' % Slepian's method which is valid everywhere disp(sprintf('\nBessel expansion method, uniformly valid\n')) % Get the Jacobi expansion coefficients [D,V,g,c]=devilliers(m,N,K,L); % What is the interval over which we are looking? defval('x',linspace(0,1,200)); if ~isnan(x) if any(x<0) error('Arguments must be zero or positive') end % Calculate the required set of Bessel functions... % ... for the normalization... [w,xGL,NGL]=gausslegendrecof(L,[],[0 1]); JGL=getscaledbessels(xGL,m,L,c); % ... for the actual points requested J=getscaledbessels(x,m,L,c); % Calculate the whole set functions as a matrix with row dimension the % evaluation points and column dimension the eigenfunction rank... % ... at the integration points EGL=JGL'*D./repmat(g,length(xGL),1); % ... at the evaluation points E=J'*D./repmat(g,length(x),1); % Prepare the normalization on the unit interval orv2=2*pi*EGL'*diag(w)*EGL; disp(sprintf('Normalization: Mean absolute error %8.3e',... mean(mean(abs(orv2-diag(diag(orv2))))))) disp(sprintf('\nGauss-Legendre normalization, Nystrom method with %i nodes',... NGL)) % Sometimes you won't want this, and undoing might be hard E2=sqrt(2*pi)*[diag(1./sqrt(diag(orv2)))*E']'; % Apply and return the normalized result E=[diag(sqrt(V(:))./sqrt(diag(orv2)))*E']'; EGL=[diag(sqrt(V(:))./sqrt(diag(orv2)))*EGL']'; else E=NaN; E2=NaN; end case 'GL' if isnan(x) error('If you only want eigenvalues select ''DV'' or ''SE''') end % By the Nystrom method on the Bessel kernel c=2*sqrt(N); % Get the weights and nodes for Gauss-Legendre integration [w,xGL,NGL]=gausslegendrecof(L,[],[0 1]); disp(sprintf('\nGauss-Legendre integration, Nystrom method with %i nodes',... NGL)) % What is the interval over which we are looking? defval('x',linspace(0,1,200)); if any(x<0) error('Arguments must be zero or positive') end switch scalem case 1 % More like Slepian, Devilliers and Shkolnisky % Evaluate the kernel at the GL points and the requested points cxxGL=c*xGL(:)*xGL(:)'; cxxint=c*x(:)*xGL(:)'; % Form the actual integration kernel KGL=sqrt(cxxGL).*besselj(m,cxxGL); Kint=sqrt(cxxint).*besselj(m,cxxint); case 0 % More like Simons, Dahlen and Wieczorek % Ahem - new stuff to be closer to our own without needing % 1/sqrt(x) - no "symmetrizing" as Slepian (1964) eq. (19) cxxGL=c*xGL(:)*xGL(:)'; cxxint=c*x(:)*xGL(:)'; KGL=sqrt(c)*besselj(m,cxxGL).*repmat(xGL(:)',length(xGL),1); Kint=sqrt(c)*besselj(m,cxxint).*repmat(xGL(:)',length(x),1); end % Get eigenfunctions of the "symmetrized" kernel - note that Matlab % may still consider them numerically asymmetric. And note that the % number of nodes appears to be a sensitive parameter. if K>=NGL-1 & K<=NGL [EGL,g]=eig(diag(sqrt(w))*KGL*diag(sqrt(w))); elseif KNGL error(sprintf('You must increase polynomial degree to at least %i',2*K-2)) end % Need to unscale E EGL=EGL./repmat(sqrt(w(:)),1,size(EGL,2)); % Now get the actual eigenvalues of the concentration problem V=diag(g).^2*c; % Sometimes you won't want this, and undoing might be hard EGL2=EGL; % Make sure the inner product over the domain is as wished EGL=EGL*sqrt(diag(V))/sqrt(2*pi); % Order actual eigenvalues and eigenfunctions downgoing [V,i]=sort(V,'descend'); V=V(:)'; EGL=EGL(:,i); EGL2=EGL2(:,i); % But usually we "interpolate" them using the same kernel at the % desired values. E=(Kint*diag(w))*EGL*inv(g); E2=(Kint*diag(w))*EGL2*inv(g); % Should build in some kind of complex/negative protection if any(imag(V)>eps) error('Complex eigenvalues - try increasing L'); else V=real(V); E=real(E); E2=real(E2); end % Check that the eigenvalues are properly contained if ~isnan(nanmean([V(V>1) V(V<0)])); V1=V-1; V1=V1(V1>0); V0=V+1; V0=V0(V0<1); if nanmean([V1 V0])>1000*eps error('Eigenvalues exceeding 0 to 1 - definitely increase L'); end end if scalem==1 || ~strcmp(method,'GL') % Check orthogonality over the interval orv2=2*pi*EGL'*diag(w)*EGL; checksit1=mean(mean(abs(orv2-diag(V)))); % Check orthonormality over the interval orv3=EGL2'*diag(w)*EGL2; checksit2=mean(mean(abs(orv3-eye(size(orv3))))); elseif scalem==0 && strcmp(method,'GL') orv2=2*pi*EGL'*diag(w)*[EGL.*repmat(xGL,1,size(EGL,2))]; % Apply and return the normalized result EGL=[diag(sqrt(V(:))./sqrt(diag(orv2)))*EGL']'; EGL2=[diag(sqrt(V(:))./sqrt(diag(orv2)))*EGL2']'; E=[diag(sqrt(V(:))./sqrt(diag(orv2)))*E']'; orv2=2*pi*EGL'*diag(w)*[EGL.*repmat(xGL,1,size(EGL,2))]; checksit1=mean(mean(abs(orv2-diag(V)))); % Check orthonormality over the interval orv3=EGL2'*diag(w)*[EGL2.*repmat(xGL,1,size(EGL2,2))]; checksit2=mean(mean(abs(orv3-eye(size(orv3))))); else error('You have run out of options, buddy') end disp(sprintf('Orthogonality of E : Mean absolute error %8.3e',... checksit1)) disp(sprintf('Orthonormality of E2: Mean absolute error %8.3e',... checksit2)) difer(checksit1,[],[],NaN) difer(checksit2,[],[],NaN) otherwise error('Specify valid method') end if ~strcmp(method,'GL') & ~isnan(x) % Check orthogonality over the interval orv2=2*pi*EGL'*diag(w)*EGL; % Sort of equal to 2*pi*nansum(E(x<=1,3).^2).*indeks(diff(x)) % Sort of equal to 2*pi*trapeze(x(x(:)<=1),[E(x<=1,3)].^2) checksit1=mean(mean(abs(orv2-diag(V)))); disp(sprintf('Orthogonality of E : Mean absolute error %8.3e',... checksit1)) difer(checksit1,[],[],NaN) end % Compare the sum over the eigenvalues with the partial Shannon number Nm=indeks(nsubm(N,m,1),m+1); disp(sprintf('Nm = %4.2g ; V1 = %4.2e',Nm,V(1))) % Only check if you've calculated at least 2N of them if K>=max(2*N,3) difer(sum(V)-Nm,8,1,'Partial Shannon numbers check out') end if ~isnan(x) % Make them start on an upswing if at all possible - might be avoided % if we wanted to keep the eigenvalue signed, which makes not much % sense of course, given the quadratic concentration problem for index=1:size(E,2) E(:,index)=E(:,index)*sign(indeks(E(~~E(:,index),index),1)); E2(:,index)=E2(:,index)*sign(indeks(E2(~~E2(:,index),index),1)); end % And remember, if you didn't like the scaling, you can undo it here if scalem==0 && ~strcmp(method,'GL') warning off MATLAB:divideByZero E=E./repmat(sqrt(x(:)),1,size(E,2)); E2=E2./repmat(sqrt(x(:)),1,size(E2,2)); disp('First element will be NaN and the FFT will be bogus') warning on MATLAB:divideByZero end end % Keep track of computation time ts=etime(clock,t0); disp(sprintf('\nComputation took %8.4f s',ts)) % Warn if things are getting tight if any(V