|
| 1 | +function [F] = chunkerflam(chnkr,kern,dval,opts) |
| 2 | +%CHUNKERFLAM build the requested FLAM compressed representation |
| 3 | +% (e.g. a recursive skeletonization factorization) of the system matrix |
| 4 | +% for given kernel and chunker description of boundary. This routine |
| 5 | +% has the same quadrature options as CHUNKERMAT. |
| 6 | +% |
| 7 | +% Syntax: sysmat = chunkerflam(chnkr,kern,opts) |
| 8 | +% |
| 9 | +% Input: |
| 10 | +% chnkr - chunker object describing boundary |
| 11 | +% kern - kernel function. By default, this should be a function handle |
| 12 | +% accepting input of the form kern(srcinfo,targinfo), where |
| 13 | +% srcinfo and targinfo are in the ptinfo struct format, i.e. |
| 14 | +% ptinfo.r - positions (2,:) array |
| 15 | +% ptinfo.d - first derivative in underlying |
| 16 | +% parameterization (2,:) |
| 17 | +% ptinfo.d2 - second derivative in underlying |
| 18 | +% parameterization (2,:) |
| 19 | +% dval - (default 0.0) float or float array. Let A be the matrix |
| 20 | +% corresponding to on-curve convolution with the provided kernel. |
| 21 | +% If a scalar is provided, the system matrix is |
| 22 | +% A + dval*eye(size(A)) |
| 23 | +% If a vector is provided, it should be length size(A,1). The |
| 24 | +% system matrix is then |
| 25 | +% A + diag(dval) |
| 26 | +% |
| 27 | +% Optional input: |
| 28 | +% opts - options structure. available options (default settings) |
| 29 | +% opts.flamtype = string ('rskelf'), type of compressed |
| 30 | +% representation to compute. Available: |
| 31 | +% |
| 32 | +% The recursive skeletonization routines |
| 33 | +% should be sufficient curves which are not |
| 34 | +% approximately space filling. |
| 35 | +% |
| 36 | +% - 'rskelf', recursive skeletonization |
| 37 | +% factorization. Can be immediately used |
| 38 | +% for system solves (via RSKELF_SV) and |
| 39 | +% determinants (via RSKELF_LOGDET) |
| 40 | +% - 'rskel', recursive skeletonization |
| 41 | +% compression. Can be immediately used for |
| 42 | +% matvecs or embedded in a sparse matrix. |
| 43 | +% Not recommended unless you have a |
| 44 | +% compelling reason. |
| 45 | +% opts.useproxy = boolean (true), use a proxy function to speed |
| 46 | +% up the FLAM compression. It may be desirable to |
| 47 | +% turn this off if there appear to be precision |
| 48 | +% issues. |
| 49 | +% |
| 50 | +% opts.quad = string ('ggqlog'), specify quadrature routine to |
| 51 | +% use. |
| 52 | +% |
| 53 | +% - 'ggqlog' uses a generalized Gaussian quadrature |
| 54 | +% designed for logarithmically singular kernels and |
| 55 | +% smooth kernels with removable singularities |
| 56 | +% - 'native' selects standard scaled Gauss-Legendre |
| 57 | +% quadrature for native functions |
| 58 | +% opts.l2scale = boolean (false), if true scale rows by |
| 59 | +% sqrt(whts) and columns by 1/sqrt(whts) |
| 60 | +% opts.occ = integer "occupancy" parameter (200) determines |
| 61 | +% how many sources are in any leaf of the tree used to |
| 62 | +% sort points. Determines the smaller dimensions of the |
| 63 | +% first set of compressions. Sent to FLAM |
| 64 | +% opts.rank_or_tol = integer or float "rank_or_tol" |
| 65 | +% parameter (1e-14). Lower precision increases speed. |
| 66 | +% Sent to FLAM |
| 67 | +% |
| 68 | +% Output: |
| 69 | +% F - the requested FLAM compressed representation of the |
| 70 | +% system matrix |
| 71 | +% |
| 72 | +% Examples: |
| 73 | +% F = chunkermat(chnkr,kern,dval); % standard options |
| 74 | +% sysmat = chunkermat(chnkr,kern,dval,opts); |
| 75 | +% |
| 76 | + |
| 77 | +F = []; |
| 78 | + |
| 79 | +if length(chnkr) > 1 |
| 80 | + chnkr = merge(chnkr); |
| 81 | +end |
| 82 | + |
| 83 | +if nargin < 3 |
| 84 | + dval = 0.0; |
| 85 | +end |
| 86 | + |
| 87 | +if nargin < 4 |
| 88 | + opts = []; |
| 89 | +end |
| 90 | + |
| 91 | +quad = 'ggqlog'; |
| 92 | +flamtype = 'rskelf'; |
| 93 | +useproxy = true; |
| 94 | +occ = 200; |
| 95 | +rank_or_tol = 1e-14; |
| 96 | + |
| 97 | + |
| 98 | +if or(chnkr.nch < 1,chnkr.k < 1) |
| 99 | + warning('empty chunker, doing nothing') |
| 100 | + sp = []; |
| 101 | + return |
| 102 | +end |
| 103 | + |
| 104 | +% determine operator dimensions using first two points |
| 105 | + |
| 106 | +srcinfo = []; targinfo = []; |
| 107 | +srcinfo.r = chnkr.r(:,1); srcinfo.d = chnkr.d(:,1); |
| 108 | +srcinfo.d2 = chnkr.d2(:,1); |
| 109 | +i2 = min(2,chnkr.npt); |
| 110 | +targinfo.r = chnkr.r(:,i2); targinfo.d = chnkr.d(:,i2); |
| 111 | +targinfo.d2 = chnkr.d2(:,i2); |
| 112 | + |
| 113 | +ftemp = kern(srcinfo,targinfo); |
| 114 | +opdims = size(ftemp); |
| 115 | + |
| 116 | +if (length(dval) == 1) |
| 117 | + dval = dval*ones(opdims(1)*chnkr.npt,1); |
| 118 | +end |
| 119 | + |
| 120 | +if (length(dval) ~= opdims(1)*chnkr.npt) |
| 121 | + warning('provided dval array is length %d. must be scalar or length %d',... |
| 122 | + length(dval),opdims(1)*chnkr.npt); |
| 123 | + return |
| 124 | +end |
| 125 | + |
| 126 | +% get opts from struct if available |
| 127 | + |
| 128 | +if isfield(opts,'quad') |
| 129 | + quad = opts.quad; |
| 130 | +end |
| 131 | +if isfield(opts,'l2scale') |
| 132 | + l2scale = opts.l2scale; |
| 133 | +end |
| 134 | + |
| 135 | +if isfield(opts,'occ') |
| 136 | + occ = opts.occ; |
| 137 | +end |
| 138 | +if isfield(opts,'rank_or_tol') |
| 139 | + rank_or_tol = opts.rank_or_tol; |
| 140 | +end |
| 141 | + |
| 142 | +% check if chosen FLAM routine implemented before doing real work |
| 143 | + |
| 144 | +if ~ (strcmpi(flamtype,'rskelf') || strcmpi(flamtype,'rskel')) |
| 145 | + warning('selected flamtype %s not available, doing nothing',flamtype) |
| 146 | + return |
| 147 | +end |
| 148 | + |
| 149 | +% get nonsmooth quadrature |
| 150 | + |
| 151 | +if strcmpi(quad,'ggqlog') |
| 152 | + |
| 153 | + type = 'log'; |
| 154 | + sp = chnk.quadggq.buildmattd(chnkr,kern,opdims,type); |
| 155 | + |
| 156 | +elseif strcmpi(quad,'native') |
| 157 | + |
| 158 | + sp = sparse(chnkr.npt,chnkr.npt); |
| 159 | + |
| 160 | +else |
| 161 | + warning('specified quadrature method not available'); |
| 162 | + return; |
| 163 | +end |
| 164 | + |
| 165 | +sp = sp + spdiags(dval,0,chnkr.npt,chnkr.npt); |
| 166 | + |
| 167 | +% prep and call flam |
| 168 | + |
| 169 | +wts = weights(chnkr); |
| 170 | +xflam = chnkr.r(:,:); |
| 171 | + |
| 172 | +width = max(max(chnkr)-min(chnkr)); |
| 173 | +optsnpxy = []; optsnpxy.rank_or_tol = rank_or_tol; |
| 174 | +optsnpxy.nsrc = occ; |
| 175 | + |
| 176 | +npxy = chnk.flam.nproxy_square(kern,width,optsnpxy); |
| 177 | + |
| 178 | +matfun = @(i,j) chnk.flam.kernbyindex(i,j,chnkr,wts,kern,opdims,sp); |
| 179 | +[pr,ptau,pw,pin] = chnk.flam.proxy_square_pts(npxy); |
| 180 | + |
| 181 | +if strcmpi(flamtype,'rskelf') |
| 182 | + ifaddtrans = true; |
| 183 | + pxyfun = @(x,slf,nbr,l,ctr) chnk.flam.proxyfun(slf,nbr,l,ctr,chnkr,wts, ... |
| 184 | + kern,opdims,pr,ptau,pw,pin,ifaddtrans); |
| 185 | + F = rskelf(matfun,xflam,occ,rank_or_tol,pxyfun); |
| 186 | +end |
| 187 | + |
| 188 | +if strcmpi(flamtype,'rskel') |
| 189 | + pxyfunr = @(rc,rx,cx,slf,nbr,l,ctr) chnk.flam.proxyfunr(rc,rx,slf,nbr,l, ... |
| 190 | + ctr,chnkr,wts,kern,opdims,pr,ptau,pw,pin); |
| 191 | + F = rskel(matfun,xflam,xflam,occ,rank_or_tol,pxyfunr); |
| 192 | +end |
| 193 | + |
| 194 | + |
| 195 | + |
| 196 | +end |
0 commit comments