Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

removes chnkr.h, closes #22 #67

Merged
merged 1 commit into from
Apr 2, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 4 additions & 5 deletions chunkie/+chnk/+quadadap/buildmat.m
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,6 @@
adj = chnkr.adj;
d = chnkr.d;
d2 = chnkr.d2;
h = chnkr.h;
n = chnkr.n;

[t,wts,u] = lege.exps(k);
Expand Down Expand Up @@ -77,7 +76,7 @@
dt = d(:,:,ibefore);
d2t = d2(:,:,ibefore);
nt = n(:,:,ibefore);
submat = chnk.adapgausswts(r,d,n,d2,h,t,bw,i,rt,dt,nt,d2t, ...
submat = chnk.adapgausswts(r,d,n,d2,t,bw,i,rt,dt,nt,d2t, ...
kern,opdims,t2,w2);

imat = 1 + (ibefore-1)*k*opdims(1);
Expand All @@ -90,7 +89,7 @@
dt = d(:,:,iafter);
nt = n(:,:,iafter);
d2t = d2(:,:,iafter);
submat = chnk.adapgausswts(r,d,n,d2,h,t,bw,i,rt,dt,nt,d2t, ...
submat = chnk.adapgausswts(r,d,n,d2,t,bw,i,rt,dt,nt,d2t, ...
kern,opdims,t2,w2);

imat = 1 + (iafter-1)*k*opdims(1);
Expand All @@ -101,7 +100,7 @@

% self

submat = chnk.quadggq.diagbuildmat(r,d,n,d2,h,[],i,kern,opdims,...
submat = chnk.quadggq.diagbuildmat(r,d,n,d2,[],i,kern,opdims,...
xs0,wts0,ainterps0kron,ainterps0);

imat = 1 + (i-1)*k*opdims(1);
Expand Down Expand Up @@ -163,7 +162,7 @@
d2t = d2(:,targfix);
nt = n(:,targfix);

submat = chnk.adapgausswts(r,d,n,d2,h,t,bw,i,rt,dt,nt,d2t, ...
submat = chnk.adapgausswts(r,d,n,d2,t,bw,i,rt,dt,nt,d2t, ...
kern,opdims,t2,w2);

imats = bsxfun(@plus,(1:opdims(1)).',opdims(1)*(targfix(:)-1).');
Expand Down
7 changes: 3 additions & 4 deletions chunkie/+chnk/+quadggq/buildmat.m
Original file line number Diff line number Diff line change
Expand Up @@ -50,7 +50,6 @@
adj = chnkr.adj;
d = chnkr.d;
d2 = chnkr.d2;
h = chnkr.h;
n = chnkr.n;

data = [];
Expand Down Expand Up @@ -101,7 +100,7 @@
if ~isempty(ilist) && ismember(ibefore,ilist) && ismember(j,ilist)
% skip construction if both chunks are in the "bad" chunk list
else
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,h,data,ibefore,j, ...
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,data,ibefore,j, ...
kern,opdims,xs1,wts1,ainterp1kron,ainterp1);

imat = 1 + (ibefore-1)*k*opdims(1);
Expand All @@ -115,7 +114,7 @@
if ~isempty(ilist) && ismember(iafter,ilist) && ismember(j,ilist)
% skip construction if both chunks are in the "bad" chunk list
else
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,h,data,iafter,j, ...
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,data,iafter,j, ...
kern,opdims,xs1,wts1,ainterp1kron,ainterp1);

imat = 1 + (iafter-1)*k*opdims(1);
Expand All @@ -129,7 +128,7 @@
if ~isempty(ilist) && ismember(j,ilist)
% skip construction if the chunk is in the "bad" chunk list
else
submat = chnk.quadggq.diagbuildmat(r,d,n,d2,h,data,j,kern,opdims,...
submat = chnk.quadggq.diagbuildmat(r,d,n,d2,data,j,kern,opdims,...
xs0,wts0,ainterps0kron,ainterps0);

imat = 1 + (j-1)*k*opdims(1);
Expand Down
7 changes: 3 additions & 4 deletions chunkie/+chnk/+quadggq/buildmattd.m
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,6 @@
adj = chnkr.adj;
d = chnkr.d;
d2 = chnkr.d2;
h = chnkr.h;
n = chnkr.n;

data = [];
Expand Down Expand Up @@ -89,7 +88,7 @@
if ~isempty(ilist) && ismember(ibefore,ilist) && ismember(j,ilist)
% skip construction if both chunks are in the "bad" chunk list
else
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,h,data,ibefore,j, ...
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,data,ibefore,j, ...
kern,opdims,xs1,wts1,ainterp1kron,ainterp1);

imat = 1 + (ibefore-1)*k*opdims(1);
Expand All @@ -106,7 +105,7 @@
if ~isempty(ilist) && ismember(iafter,ilist) && ismember(j,ilist)
% skip construction if both chunks are in the "bad" chunk list
else
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,h,data,iafter,j, ...
submat = chnk.quadggq.nearbuildmat(r,d,n,d2,data,iafter,j, ...
kern,opdims,xs1,wts1,ainterp1kron,ainterp1);


Expand All @@ -124,7 +123,7 @@
if ~isempty(ilist) && ismember(j,ilist)
% skip construction if the chunk is in the "bad" chunk list
else
submat = chnk.quadggq.diagbuildmat(r,d,n,d2,h,data,j,kern,opdims,...
submat = chnk.quadggq.diagbuildmat(r,d,n,d2,data,j,kern,opdims,...
xs0,wts0,ainterps0kron,ainterps0);

imat = 1 + (j-1)*k*opdims(1);
Expand Down
6 changes: 3 additions & 3 deletions chunkie/+chnk/+quadggq/diagbuildmat.m
Original file line number Diff line number Diff line change
@@ -1,10 +1,10 @@
function submat = diagbuildmat(r,d,n,d2,h,data,i,fkern,opdims,...
function submat = diagbuildmat(r,d,n,d2,data,i,fkern,opdims,...
xs0,whts0,ainterps0kron,ainterps0)
%CHNK.QUADGGQ.DIAGBUILDMAT

% grab specific boundary data

rs = r(:,:,i); ds = d(:,:,i); d2s = d2(:,:,i); hs = h(i);
rs = r(:,:,i); ds = d(:,:,i); d2s = d2(:,:,i);
ns = n(:,:,i);
if(isempty(data))
dd = [];
Expand All @@ -29,7 +29,7 @@
d2fine{i} = (ainterps0{i}*(d2s.')).';
dfinenrm = sqrt(sum(dfine{i}.^2,1));
nfine{i} = [dfine{i}(2,:); -dfine{i}(1,:)]./dfinenrm;
dsdt{i} = (dfinenrm(:)).*whts0{i}*hs;
dsdt{i} = (dfinenrm(:)).*whts0{i};
end

srcinfo = [];
Expand Down
6 changes: 3 additions & 3 deletions chunkie/+chnk/+quadggq/nearbuildmat.m
Original file line number Diff line number Diff line change
@@ -1,12 +1,12 @@
function submat = nearbuildmat(r,d,n,d2,h,data,i,j,fkern,opdims,...
function submat = nearbuildmat(r,d,n,d2,data,i,j,fkern,opdims,...
xs1,whts1,ainterp1kron,ainterp1)
%CHNKR.QUADGGQ.NEARBUILDMAT

% grab specific boundary data

rs = r(:,:,j); ds = d(:,:,j); d2s = d2(:,:,j); ns = n(:,:,j);
rt = r(:,:,i); dt = d(:,:,i); d2t = d2(:,:,i); nt = n(:,:,i);
hs = h(j);

if(isempty(data))
dd = [];
else
Expand Down Expand Up @@ -64,7 +64,7 @@

dfinenrm = sqrt(sum(dfine.^2,1));
%dfinenrm = dfine(1,:,:); % for complex contour, by SJ 09/30/21
dsdt = dfinenrm(:).*whts1(:)*hs;
dsdt = dfinenrm(:).*whts1(:);

dsdtndim2 = repmat(dsdt(:).',opdims(2),1);
dsdtndim2 = dsdtndim2(:);
Expand Down
5 changes: 2 additions & 3 deletions chunkie/+chnk/+quadnative/buildmat.m
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,6 @@
r = chnkr.rstor;
d = chnkr.dstor;
d2 = chnkr.d2stor;
h = chnkr.hstor;
n = chnkr.nstor;

[dim,k,~] = size(r);
Expand All @@ -34,10 +33,10 @@

srcinfo = []; srcinfo.r = rs; srcinfo.d = ds; srcinfo.d2 = d2s; srcinfo.n = ns;
targinfo = []; targinfo.r = rt; targinfo.d = dt; targinfo.d2 = d2t; targinfo.n = nt;
hs = h(j); ht = h(i);


dsnrms = sqrt(sum(ds.^2,1));
ws = kron(hs(:),wts(:));
ws = repmat(wts(:),length(j),1);

dsdt = dsnrms(:).*ws;

Expand Down
68 changes: 0 additions & 68 deletions chunkie/+chnk/+rcip/Rcomp.m

This file was deleted.

11 changes: 5 additions & 6 deletions chunkie/+chnk/+rcip/Rcompchunk.m
Original file line number Diff line number Diff line change
Expand Up @@ -106,7 +106,6 @@
d2 = chnkri.d2(:,:,ie);
il = chnkri.adj(1,ie);
ir = chnkri.adj(2,ie);
h = chnkri.h(ie);
if (il > 0 && ir < 0)
nextchunk(i) = il;
ileftright(i) = 1;
Expand All @@ -116,8 +115,8 @@
rcs(:,:,i) = sbcrmat*(r.');
dcs(:,:,i) = u*(d.');
d2cs(:,:,i) = u*(d2.');
dscal(i) = h*2;
d2scal(i) = h^2*4;
dscal(i) = 2;
d2scal(i) = 4;
elseif (il < 0 && ir > 0)
nextchunk(i) = ir;
ileftright(i) = -1;
Expand All @@ -126,8 +125,8 @@
rcs(:,:,i) = sbclmat*(r.');
dcs(:,:,i) = u*(d.');
d2cs(:,:,i) = u*(d2.');
dscal(i) = h*2;
d2scal(i) = h^2*4;
dscal(i) = 2;
d2scal(i) = 4;
else
error('RCIP: edge chunk not adjacent to one vertex and one neighbor')
end
Expand Down Expand Up @@ -197,7 +196,7 @@
chnkrlocal(i).r(:,:,nchi+1) = chnkr(ic).r(:,:,nc)-ctr(:,i);
chnkrlocal(i).d(:,:,nchi+1) = chnkr(ic).d(:,:,nc);
chnkrlocal(i).d2(:,:,nchi+1) = chnkr(ic).d2(:,:,nc);
chnkrlocal(i).h(nchi+1) = chnkr(ic).h(nc);

if ileftright(i) == -1
chnkrlocal(i).adj(1,nchi+1) = nchi;
chnkrlocal(i).adj(2,nchi+1) = -1;
Expand Down
7 changes: 4 additions & 3 deletions chunkie/+chnk/+rcip/chunkerfunclocal.m
Original file line number Diff line number Diff line change
Expand Up @@ -71,10 +71,11 @@

ts = a + (b-a)*(xs+1)/2;
[out{:}] = fcurve(ts);
h = (b-a)/2;
chnkr.rstor(:,:,i) = reshape(out{1},dim,k);
chnkr.dstor(:,:,i) = reshape(out{2},dim,k);
chnkr.d2stor(:,:,i) = reshape(out{3},dim,k);
chnkr.hstor(i) = (b-a)/2;
chnkr.dstor(:,:,i) = reshape(out{2},dim,k)*h;
chnkr.d2stor(:,:,i) = reshape(out{3},dim,k)*h*h;

end

chnkr.adjstor(:,1:nch) = adjs(:,1:nch);
Expand Down
6 changes: 1 addition & 5 deletions chunkie/+chnk/adapgausskerneval.m
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
function [fints,maxrecs,numints,iers] = adapgausskerneval(r,d,n,d2,h,ct,bw,j,...
function [fints,maxrecs,numints,iers] = adapgausskerneval(r,d,n,d2,ct,bw,j,...
dens,rt,nt,dt,d2t,kern,opdims,t,w,opts)
%CHNK.ADAPGAUSSKERNEVAL adaptive integration for interaction of kernel on chunk
% at targets
Expand All @@ -13,7 +13,6 @@
% r - chnkr nodes
% d - chnkr derivatives at nodes
% d2 - chnkr 2nd derivatives at nodes
% h - lengths of chunks in parameter space
% ct - Legendre nodes at order of chunker
% bw - barycentric interpolation weights for Legendre nodes at order of
% chunker
Expand Down Expand Up @@ -66,7 +65,6 @@
ds = d(:,:,j);
ns = n(:,:,j);
d2s = d2(:,:,j);
hs = h(j);
jstart = opdims(2)*k*(j-1)+1;
jend = opdims(2)*k*j;
densj = reshape(dens(jstart:jend),opdims(2),k);
Expand Down Expand Up @@ -161,8 +159,6 @@

end

fints = fints*hs;

end

function val = oneintp(a,b,rs,ds,ns,d2s,densj,ct,bw, ...
Expand Down
5 changes: 1 addition & 4 deletions chunkie/+chnk/adapgausswts.m
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
function [mat,maxrecs,numints,iers] = adapgausswts(r,d,n,d2,h,ct,bw,j,...
function [mat,maxrecs,numints,iers] = adapgausswts(r,d,n,d2,ct,bw,j,...
rt,dt,nt,d2t,kern,opdims,t,w,opts)
%CHNK.ADAPGAUSSWTS adaptive integration for interaction of kernel on chunk
% at targets
Expand Down Expand Up @@ -68,7 +68,6 @@
ds = d(:,:,j);
ns = n(:,:,j);
d2s = d2(:,:,j);
hs = h(j);

stack = zeros(2,maxdepth);
vals = zeros(opdims(1)*opdims(2)*k,maxdepth);
Expand Down Expand Up @@ -167,8 +166,6 @@

end

mat = mat*hs;

end

function val = oneintp(a,b,rs,ds,ns,d2s,ct,bw,rt,dt,nt,d2t,kern,opdims,t,w)
Expand Down
9 changes: 4 additions & 5 deletions chunkie/+chnk/pquadwts.m
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
function mat = pquadwts(r,d,n,d2,h,ct,bw,j,...
function mat = pquadwts(r,d,n,d2,ct,bw,j,...
rt,dt,nt,d2t,kern,opdims,t,w,opts,intp_ab,intp)
%CHNK.INTERPQUADWTS product integration for interaction of kernel on chunk
% at targets
Expand All @@ -14,7 +14,6 @@
% d - chnkr derivatives at nodes
% n - chnkr normals at nodes
% d2 - chnkr 2nd derivatives at nodes
% h - lengths of chunks in parameter space
% ct - Legendre nodes at order of chunker
% bw - barycentric interpolation weights for Legendre nodes at order of
% chunker
Expand All @@ -34,11 +33,11 @@
%

% Helsing-Ojala (interior/exterior?)
h_i = h(j);

xlohi = intp_ab*(r(1,:,j)'+1i*r(2,:,j)'); % panel end points
r_i = intp*(r(1,:,j)'+1i*r(2,:,j)'); % upsampled panel
d_i = h_i*(intp*(d(1,:,j)'+1i*d(2,:,j)')); % r'
d2_i = h_i^2*(intp*(d2(1,:,j)'+1i*d2(2,:,j)')); % r''
d_i = (intp*(d(1,:,j)'+1i*d(2,:,j)')); % r'
d2_i = (intp*(d2(1,:,j)'+1i*d2(2,:,j)')); % r''
sp = abs(d_i); tang = d_i./sp; % speed, tangent
n_i = -1i*tang; % normal
cur = -real(conj(d2_i).*n_i)./sp.^2; % curvature
Expand Down
Loading
Loading