source: trunk/src/civ_uvmat.m @ 243

Last change on this file since 243 was 243, checked in by sommeria, 13 years ago

bug repair in uvmat

File size: 19.8 KB
Line 
1% To develop....
2function [Data,errormsg]= civ_uvmat(Param,ncfile)
3Data.ListGlobalAttribute={'Conventions','Program','CivStage'};
4Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
5Data.Program='civ_uvmat';
6Data.CivStage=0;%default
7
8%% Civ1
9if isfield (Param,'Civ1')
10    par_civ1=Param.Civ1;
11    str2num(par_civ1.rho)
12    image1=imread(par_civ1.filename_ima_a);
13    image2=imread(par_civ1.filename_ima_b);
14    stepx=str2num(par_civ1.dx);
15    stepy=str2num(par_civ1.dy);
16    ibx2=ceil(str2num(par_civ1.ibx)/2);
17    iby2=ceil(str2num(par_civ1.iby)/2);
18    isx2=ceil(str2num(par_civ1.isx)/2);
19    isy2=ceil(str2num(par_civ1.isy)/2);
20    shiftx=str2num(par_civ1.shiftx);
21    shifty=str2num(par_civ1.shifty);
22    miniy=max(1+isy2-shifty,1+iby2);
23    minix=max(1+isx2-shiftx,1+ibx2);
24    maxiy=min(size(image1,1)-isy2-shifty,size(image1,1)-iby2);
25    maxix=min(size(image1,2)-isx2-shiftx,size(image1,2)-ibx2);
26    [GridX,GridY]=meshgrid(minix:stepx:maxix,miniy:stepy:maxiy);
27    PointCoord(:,1)=reshape(GridX,[],1);
28    PointCoord(:,2)=reshape(GridY,[],1);
29   
30    % caluclate velocity data (y and v in indices, reverse to y component)
31    [xtable ytable utable vtable ctable F] = pivlab (image1,image2,ibx2,iby2,isx2,isy2,shiftx,shifty,PointCoord,str2num(par_civ1.rho), []);
32    list_param=(fieldnames(par_civ1))';
33    list_remove={'pxcmx','pxcmy','npx','npy','gridflag','maskflag','term_a','term_b','T0'};
34    index=zeros(size(list_param));
35    for ilist=1:length(list_remove)
36        index=strcmp(list_remove{ilist},list_param);
37        if ~isempty(find(index,1))
38            list_param(index)=[];
39        end
40    end
41    for ilist=1:length(list_param)
42        Civ1_param{ilist}=['Civ1_' list_param{ilist}];
43        eval(['Data.Civ1_' list_param{ilist} '=Param.Civ1.' list_param{ilist} ';'])
44    end
45    if isfield(Data,'Civ1_gridname') && strcmp(Data.Civ1_gridname(1:6),'noFile')
46        Data.Civ1_gridname='';
47    end
48    if isfield(Data,'Civ1_maskname') && strcmp(Data.Civ1_maskname(1:6),'noFile')
49        Data.Civ1_maskname='';
50    end
51    Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ1_param {'Civ1_Time','Civ1_Dt'}];
52    Data.Civ1_Time=str2double(par_civ1.T0);
53    Data.Civ1_Dt=str2double(par_civ1.Dt);
54    Data.ListVarName={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F'};%  cell array containing the names of the fields to record
55    Data.VarDimName={'nbvec1','nbvec1','nbvec1','nbvec1','nbvec1','nbvec1'};
56    Data.VarAttribute{1}.Role='coord_x';
57    Data.VarAttribute{2}.Role='coord_y';
58    Data.VarAttribute{3}.Role='vector_x';
59    Data.VarAttribute{4}.Role='vector_y';
60    Data.VarAttribute{5}.Role='warnflag';
61    Data.Civ1_X=reshape(xtable,[],1);
62    Data.Civ1_Y=reshape(size(image1,1)-ytable+1,[],1);
63    Data.Civ1_U=reshape(utable,[],1);
64    Data.Civ1_V=reshape(-vtable,[],1);
65    Data.Civ1_C=reshape(ctable,[],1);
66    Data.Civ1_F=reshape(F,[],1);
67    Data.CivStage=1;
68else
69    Data=nc2struct(ncfile)%read existing netcdf file
70    if isfield(Data,'absolut_time_T0')%read civx file
71        var={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F','Civ1_FF';'vec_X','vec_Y','vec_U','vec_V','vec_C','vec_F','vec_FixFlag'};
72        %var=varcivx_generator('velocity','Civ1');%determine the names of constants and variables to read
73        [Data,vardetect,ichoice]=nc2struct(ncfile,var);%read the variables in the netcdf file
74        Data.ListGlobalAttribute=[{'Conventions','Program','CivStage'} Data.ListGlobalAttribute {'Civ1_Time','Civ1_Dt'}];
75        Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
76        Data.Program='civ_uvmat';
77        Data.Civ1_Time=double(Data.absolut_time_T0);
78        Data.Civ1_Dt=double(Data.dt);
79        Data.VarDimName={'nbvec1','nbvec1','nbvec1','nbvec1','nbvec1','nbvec1','nbvec1'};
80        Data.VarAttribute{1}.Role='coord_x';
81        Data.VarAttribute{2}.Role='coord_y';
82        Data.VarAttribute{3}.Role='vector_x';
83        Data.VarAttribute{4}.Role='vector_y';
84        Data.VarAttribute{5}.Role='ancillary';
85        Data.VarAttribute{6}.Role='warnflag';
86        Data.VarAttribute{7}.Role='errorflag';
87        Data.CivStage=1;
88    end
89end
90
91%% Fix1
92if isfield (Param,'Fix1')
93    ListFixParam=fieldnames(Param.Fix1);
94    for ilist=1:length(ListFixParam)
95        ParamName=ListFixParam{ilist};
96        ListName=['Fix1_' ParamName];
97        ['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];']
98        eval(['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];'])
99        eval(['Data.' ListName '=Param.Fix1.' ParamName ';'])
100    end
101%     Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Fix1_WarnFlags','Fix1_TreshCorr','Fix1_TreshVel','Fix1_UpperBoundTest'}];
102%     Data.Fix1_WarnFlags=Param.Fix1.WarnFlags;
103%     Data.Fix1_ThreshCorr=Param.Fix1.ThreshCorr;
104%     Data.Fix1_ThreshVel=Param.Fix1.ThreshVel;
105%     Data.Fix1_UpperBoundTest=Param.Fix1.UpperBoundTest;
106    Data.ListVarName=[Data.ListVarName {'Civ1_FF'}];
107    Data.VarDimName=[Data.VarDimName {'nbvec1'}];
108    nbvar=length(Data.ListVarName);
109    Data.VarAttribute{nbvar}.Role='errorflag';
110    [Data.Civ1_FF]=fix_uvmat(Param.Fix1,Data.Civ1_F,Data.Civ1_C,Data.Civ1_U,Data.Civ1_V);
111    Data.CivStage=2;                               
112end   
113%% Patch1
114if isfield (Param,'Patch1')
115    Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Patch1_Rho','Patch1_Threshold','Patch1_SubDomain'}];
116    Data.Patch1_Rho=str2double(Param.Patch1.Rho);
117    Data.Patch1_Threshold=str2double(Param.Patch1.Threshold);
118    Data.Patch1_SubDomain=str2double(Param.Patch1.SubDomain);
119    Data.ListVarName=[Data.ListVarName {'Patch1_U','Patch1_V'}];
120    Data.VarDimName=[Data.VarDimName {'nbvec1','nbvec1'}];
121    nbvar=length(Data.ListVarName);
122    Data.VarAttribute{nbvar-1}.Role='vector_x';
123    Data.VarAttribute{nbvar}.Role='vector_y';
124    Data.Patch1_U=zeros(size(Data.Civ1_X));
125    Data.Patch1_V=zeros(size(Data.Civ1_X));
126    if isfield(Data,'Civ1_FF')
127        ind_good=find(Data.Civ1_FF==0);
128    else
129        ind_good=1:numel(Data.Civ1_X);
130    end
131    Data.Civ1_X
132    [Ures, Vres,FFres]=...
133                            patch_uvmat(Data.Civ1_X(ind_good)',Data.Civ1_Y(ind_good)',Data.Civ1_U(ind_good)',Data.Civ1_V(ind_good)',Data.Patch1_Rho,Data.Patch1_Threshold,Data.Patch1_SubDomain);
134                        size(Ures)
135                        size(Vres)
136                        size(FFres)
137                        size(ind_good)
138      Data.Patch1_U(ind_good)=Ures;
139      Data.Patch1_V(ind_good)=Vres;
140      Data.Civ1_FF(ind_good)=FFres
141      Data.CivStage=3;                               
142end   
143%% write result
144errormsg=struct2nc(ncfile,Data);
145
146
147
148%'RUN_FIX': function for fixing velocity fields:
149%-----------------------------------------------
150% RUN_FIX(filename,field,flagindex,thresh_vecC,thresh_vel,iter,flag_mask,maskname,fileref,fieldref)
151%
152%filename: name of the netcdf file (used as input and output)
153%field: structure specifying the names of the fields to fix (depending on civ1 or civ2)
154    %.vel_type='civ1' or 'civ2';
155    %.nb=name of the dimension common to the field to fix ('nb_vectors' for civ1);
156    %.fixflag=name of fix flag variable ('vec_FixFlag' for civ1)
157%flagindex: flag specifying which values of vec_f are removed:
158        % if flagindex(1)=1: vec_f=-2 vectors are removed
159        % if flagindex(2)=1: vec_f=3 vectors are removed
160        % if flagindex(3)=1: vec_f=2 vectors are removed (if iter=1) or vec_f=4 vectors are removed (if iter=2)
161%iter=1 for civ1 fields and iter=2 for civ2 fields
162%thresh_vecC: threshold in the image correlation vec_C
163%flag_mask: =1 mask used to remove vectors (0 else)
164%maskname: name of the mask image file for fix
165%thresh_vel: threshold on velocity, or on the difference with the reference file fileref if exists
166%inf_sup=1: remove values smaller than threshold thresh_vel, =2, larger than threshold
167%fileref: .nc file name for a reference velocity (='': refrence 0 used)
168%fieldref: 'civ1','filter1'...feld used in fileref
169
170function FF=fix_uvmat(Param,F,C,U,V)
171%error=[]; %default
172FF=zeros(size(F));%default
173
174%criterium on warn flags
175if isfield (Param,'WarnFlags')
176    for iflag=1:numel(Param.WarnFlags)
177        FF=(FF==1| F==Param.WarnFlags(iflag));
178    end
179end
180
181%criterium on correlation values
182if isfield (Param,'LowerBoundCorr')
183    FF=FF==1 | C<Param.LowerBoundCorr;
184end
185
186if isfield (Param,'LowerBoundVel')&& ~isequal(Param.LowerBoundVel,0)
187    thresh=Param.LowerBoundVel*Param.LowerBoundVel;
188    FF=FF==1 | (U.*U+V.*V)<thresh;
189end
190if isfield (Param,'UpperBoundVel')&& ~isequal(Param.UpperBoundVel,0)
191    thresh=Param.UpperBoundVel*Param.UpperBoundVel;
192    FF=FF==1 | (U.*U+V.*V)>thresh;
193end
194 
195%
196% % criterium on velocity values
197% delta_u=Field.U;%default without ref file
198% delta_v=Field.V;
199% if exist('fileref','var') && ~isempty(fileref)
200%     if ~exist(fileref,'file')
201%         error='reference file not found in RUN_FIX.m';
202%         display(error);
203%         return
204%     end
205%     FieldRef=read_civxdata(fileref,[],fieldref);   
206%     if isfield(FieldRef,'FF')
207%         index_true=find(FieldRef.FF==0);
208%         FieldRef.X=FieldRef.X(index_true);
209%         FieldRef.Y=FieldRef.Y(index_true);
210%         FieldRef.U=FieldRef.U(index_true);
211%         FieldRef.V=FieldRef.V(index_true);
212%     end
213%     if ~isfield(FieldRef,'X') || ~isfield(FieldRef,'Y') || ~isfield(FieldRef,'U') || ~isfield(FieldRef,'V')
214%         error='reference file is not a velocity field in RUN_FIX.m '; %bad input file
215%         return
216%     end
217%     if length(FieldRef.X)<=1
218%         errordlg('reference field with one vector or less in RUN_FIX.m')
219%         return
220%     end
221%     vec_U_ref=griddata_uvmat(FieldRef.X,FieldRef.Y,FieldRef.U,Field.X,Field.Y);  %interpolate vectors in the ref field
222%     vec_V_ref=griddata_uvmat(FieldRef.X,FieldRef.Y,FieldRef.V,Field.X,Field.Y);  %interpolate vectors in the ref field to the positions  of the main field     
223%     delta_u=Field.U-vec_U_ref;%take the difference with the interpolated ref field
224%     delta_v=Field.V-vec_V_ref;
225% end
226% thresh_vel_x=thresh_vel;
227% thresh_vel_y=thresh_vel;
228% if isequal(inf_sup,1)
229%     flag5=abs(delta_u)<thresh_vel_x & abs(delta_v)<thresh_vel_y &(flag1~=1)&(flag2~=1)&(flag3~=1)&(flag4~=1);
230% elseif isequal(inf_sup,2)
231%     flag5=(abs(delta_u)>thresh_vel_x | abs(delta_v)>thresh_vel_y) &(flag1~=1)&(flag2~=1)&(flag3~=1)&(flag4~=1);
232% end
233%
234%             % flag7 introduce a grey mask, matrix M
235% if isequal (flag_mask,1)
236%    M=imread(maskname);
237%    nxy=size(M);
238%    M=reshape(M,1,nxy(1)*nxy(2));
239%    rangx0=[0.5 nxy(2)-0.5];
240%    rangy0=[0.5 nxy(1)-0.5];
241%    vec_x1=Field.X-Field.U/2;%beginning points
242%    vec_x2=Field.X+Field.U/2;%end points of vectors
243%    vec_y1=Field.Y-Field.V/2;%beginning points
244%    vec_y2=Field.Y+Field.V/2;%end points of vectors
245%    indx=1+round((nxy(2)-1)*(vec_x1-rangx0(1))/(rangx0(2)-rangx0(1)));% image index x at abcissa vec_x
246%    indy=1+round((nxy(1)-1)*(vec_y1-rangy0(1))/(rangy0(2)-rangy0(1)));% image index y at ordinate vec_y   
247%    test_in=~(indx < 1 |indy < 1 | indx > nxy(2) |indy > nxy(1)); %=0 out of the mask image, 1 inside
248%    indx=indx.*test_in+(1-test_in); %replace indx by 1 out of the mask range
249%    indy=indy.*test_in+(1-test_in); %replace indy by 1 out of the mask range
250%    ICOMB=((indx-1)*nxy(1)+(nxy(1)+1-indy));%determine the indices in the image reshaped in a Matlab vector
251%    Mvalues=M(ICOMB);
252%    flag7b=((20 < Mvalues) & (Mvalues < 200))| ~test_in';
253%    indx=1+round((nxy(2)-1)*(vec_x2-rangx0(1))/(rangx0(2)-rangx0(1)));% image index x at abcissa Field.X
254%    indy=1+round((nxy(1)-1)*(vec_y2-rangy0(1))/(rangy0(2)-rangy0(1)));% image index y at ordinate vec_y
255%    test_in=~(indx < 1 |indy < 1 | indx > nxy(2) |indy > nxy(1)); %=0 out of the mask image, 1 inside
256%    indx=indx.*test_in+(1-test_in); %replace indx by 1 out of the mask range
257%    indy=indy.*test_in+(1-test_in); %replace indy by 1 out of the mask range
258%    ICOMB=((indx-1)*nxy(1)+(nxy(1)+1-indy));%determine the indices in the image reshaped in a Matlab vector
259%    Mvalues=M(ICOMB);
260%    flag7e=((Mvalues > 20) & (Mvalues < 200))| ~test_in';
261%    flag7=(flag7b|flag7e)';
262% else
263%    flag7=0;
264% end   
265% flagmagenta=flag1|flag2|flag3|flag4|flag5|flag7;
266% fixflag_unit=Field.FF-10*floor(Field.FF/10); %unity term of fix_flag
267
268
269
270%------------------------------------------------------------------------
271% patch function
272function [U_patch,V_patch,FF,SubRangx,SubRangy,X_ctrs,Y_ctrs] =patch_uvmat(X,Y,U,V,Rho,Threshold,SubDomain)
273%subdomain decomposition
274warning off
275U=reshape(U,[],1);
276V=reshape(V,[],1);
277X=reshape(X,[],1);
278Y=reshape(Y,[],1);
279nbvec=numel(X);
280NbSubDomain=ceil(nbvec/SubDomain)
281MinX=min(X)
282MinY=min(Y)
283MaxX=max(X)
284MaxY=max(Y)
285RangX=MaxX-MinX;
286RangY=MaxY-MinY;
287AspectRatio=RangY/RangX
288NbSubDomainX=ceil(sqrt(NbSubDomain/AspectRatio))
289NbSubDomainY=ceil(sqrt(NbSubDomain*AspectRatio))
290
291SizX=RangX/NbSubDomainX;%width of subdomains
292SizY=RangY/NbSubDomainY;%height of subdomains
293CentreX=linspace(MinX+SizX/2,MaxX-SizX/2,NbSubDomainX);
294CentreY=linspace(MinY+SizY/2,MaxY-SizY/2,NbSubDomainY);
295SubIndexX=ceil((X-MinX)/SizX);%subdomain index of vectors
296SubIndexY=ceil((Y-MinY)/SizY);
297rho=RangX*RangY*Rho;%optimum rho increase as teh area of the subdomain
298U_tps=zeros(length(X),NbSubDomainY,NbSubDomainX);%default spline
299V_tps=zeros(length(X),NbSubDomainY,NbSubDomainX);%default spline
300X_dist=zeros(length(X),NbSubDomainY,NbSubDomainX);%default spline
301Y_dist=zeros(length(X),NbSubDomainY,NbSubDomainX);%default spline
302U_smooth=zeros(length(X),NbSubDomainY,NbSubDomainX);
303V_smooth=zeros(length(X),NbSubDomainY,NbSubDomainX);
304dist_ctre=zeros(length(X),NbSubDomainY,NbSubDomainX);
305FF=zeros(length(X),1);
306for isubx=1:NbSubDomainX
307    for isuby=1:NbSubDomainY
308        SubRangx(isuby,isubx,:)=[CentreX(isubx)-SizX/2 CentreX(isubx)+SizX/2];
309        SubRangy(isuby,isubx,:)=[CentreY(isubx)-SizY/2 CentreY(isubx)+SizY/2];
310        for iter=1:3 %increase the subdomain during three iterations at most
311            ind_sel=find(X>SubRangx(isuby,isubx,1) & X<SubRangx(isuby,isubx,2) & Y>SubRangy(isuby,isubx,1) & Y<SubRangy(isuby,isubx,2)); 
312            size(ind_sel)
313            if numel(ind_sel)<SubDomain/4;% too few selected vectors, increase the subrange for next iteration
314                SubRangx(isuby,isubx,1)=SubRangx(isuby,isubx,1)-SizX/4;
315                SubRangx(isuby,isubx,2)=SubRangx(isuby,isubx,2)+SizX/4;
316                SubRangy(isuby,isubx,1)=SubRangy(isuby,isubx,1)-SizY/4;
317                SubRangy(isuby,isubx,2)=SubRangy(isuby,isubx,2)+SizY/4;
318            else
319                [U_smooth_sub,U_tps_sub]=tps_uvmat(X(ind_sel),Y(ind_sel),U(ind_sel),rho);
320                [V_smooth_sub,V_tps_sub]=tps_uvmat(X(ind_sel),Y(ind_sel),V(ind_sel),rho);
321                size(U_smooth_sub)
322                size(U(ind_sel))
323                UDiff=U_smooth_sub-U(ind_sel);
324                VDiff=V_smooth_sub-V(ind_sel);
325                NormDiff=UDiff.*UDiff+VDiff.*VDiff;
326                FF(ind_sel)=20*(NormDiff>Threshold);%put FF value to 20 to identify the criterium of elimmination
327                ind_ind_sel=find(FF(ind_sel)==0);
328                if isequal(numel(ind_ind_sel),numel(ind_sel))
329                    U_smooth(ind_sel,isuby,isubx)=U_smooth_sub;
330                    V_smooth(ind_sel,isuby,isubx)=V_smooth_sub;
331                    break
332                elseif numel(ind_ind_sel)<SubDomain/4;% too few selected vectors, increase the subrange for next iteration
333                    SubRangx(isuby,isubx,1)=SubRangx(isuby,isubx,1)-SizX/4;
334                    SubRangx(isuby,isubx,2)=SubRangx(isuby,isubx,2)+SizX/4;
335                    SubRangy(isuby,isubx,1)=SubRangy(isuby,isubx,1)-SizY/4;
336                    SubRangy(isuby,isubx,2)=SubRangy(isuby,isubx,2)+SizY/4;
337                else
338                    [U_smooth(ind_sel(ind_ind_sel),isuby,isubx),U_tps(ind_sel(ind_ind_sel),isuby,isubx),U_tps3(isuby,isubx)]=tps_uvmat(X(ind_sel(ind_ind_sel)),Y(ind_sel(ind_ind_sel)),U(ind_sel(ind_ind_sel)),rho);
339                    [V_smooth(ind_sel(ind_ind_sel),isuby,isubx),V_tps(ind_sel(ind_ind_sel),isuby,isubx),V_tps3(isuby,isubx)]=tps_uvmat(X(ind_sel(ind_ind_sel)),Y(ind_sel(ind_ind_sel)),V(ind_sel(ind_ind_sel)),rho);
340                     break                 
341                end
342            end       
343        end       
344        X_ctrs(isuby,isubx)=mean(X(ind_sel));%gravity centre of selected points
345        Y_ctrs(isuby,isubx)=mean(Y(ind_sel));%positions of tps sources for the subdomain i,     
346        dist_ctre(ind_sel,isuby,isubx)=sqrt(abs(((X(ind_sel)-X_ctrs(isuby,isubx)).*(X(ind_sel)-X_ctrs(isuby,isubx)))+((Y(ind_sel)-Y_ctrs(isuby,isubx)).*(Y(ind_sel)-Y_ctrs(isuby,isubx)))));
347    end
348end
349U_patch=sum(sum(U_smooth.*dist_ctre,3),2)./sum(sum(dist_ctre,3),2);
350V_patch=sum(sum(V_smooth.*dist_ctre,3),2)./sum(sum(dist_ctre,3),2);
351
352%------------------------------------------------------------------------
353%fasshauer@iit.edu MATH 590 ? Chapter 19 32
354% X,Y initial coordiantes
355% XI vector, YI column vector for the grid of interpolation points
356function [U_smooth,U_tps,U_tps3]=tps_uvmat(X,Y,U,rho)
357%------------------------------------------------------------------------
358%rho smoothing parameter
359ep = 1;
360X=reshape(X,[],1);
361Y=reshape(Y,[],1);
362rhs = reshape(U,[],1);
363% if exist('FF','var')
364% test_false=isnan(rhs)|FF~=0;
365% else
366%     test_false=isnan(rhs);
367% end
368% X(test_false)=[];
369% Y(test_false)=[];
370% rhs(test_false)=[];
371%randn('state',3); rhs = rhs + 0.03*randn(size(rhs));
372rhs = [rhs; zeros(3,1)];
373dsites = [X Y];% coordinates of measurement sites
374ctrs = dsites;%radial base functions are located at the measurement sites
375DM_data = DistanceMatrix(dsites,ctrs);%2D matrix of distances between spline centres (=initial points) ctrs
376% if size(XI,1)==1 && size(YI,2)==1 % XI vector, YI column vector
377%      [XI,YI]=meshgrid(XI,YI);
378% end
379% [npy,npx]=size(XI);
380% epoints = [reshape(XI,[],1) reshape(YI,[],1)];
381IM_sites = tps(ep,DM_data);%values of thin plate at site points
382IM = IM_sites + rho*eye(size(IM_sites));%  rho=1/(2*omega) , omega given by fasshauer;
383PM=[ones(size(dsites,1),1) dsites];
384IM=[IM PM; [PM' zeros(3,3)]];
385%fprintf('Condition number estimate: %e\n',condest(IM))
386%DM_eval = DistanceMatrix(epoints,ctrs);%2D matrix of distances between extrapolation points epoints and spline centres (=site points) ctrs
387%EM = tps(ep,DM_eval);%values of thin plate
388%PM = [ones(size(epoints,1),1) epoints];
389%EM = [EM PM];
390U_tps=(IM\rhs);
391PM = [ones(size(dsites,1),1) dsites];
392EM = [IM_sites PM];
393U_smooth=EM *U_tps;
394U_tps3=U_tps(end-2:end);
395U_tps=U_tps(1:end-3);
396
397
398% U_patch = EM * spline_coeff;
399% U_patch=reshape(U_patch,npy,npx);
400% PM = [ones(size(dsites,1),1) dsites];
401% EM = [IM_sites PM];
402% U(test_false)=[];
403% U_nodes=EM * spline_coeff;
404
405%exact = testfunctions(epoints);
406%maxerr = norm(Pf-exact,inf);
407% PlotSurf(xe,ye,Pf,neval,exact,maxerr,[160,20]);
408% PlotError2D(xe,ye,Pf,exact,maxerr,neval,[160,20]);
409
410
411
412  % DM = DistanceMatrix(dsites,ctrs)
413% Forms the distance matrix of two sets of points in R^s,
414% i.e., DM(i,j) = || datasite_i - center_j ||_2.
415% Input
416%   dsites: Mxs matrix representing a set of M data sites in R^s
417%              (i.e., each row contains one s-dimensional point)
418%   ctrs:   Nxs matrix representing a set of N centers in R^s
419%              (one center per row)
420% Output
421
422
423%   DM:     MxN matrix whose i,j position contains the Euclidean
424%              distance between the i-th data site and j-th center
425  function DM = DistanceMatrix(dsites,ctrs)
426  [M,s] = size(dsites); [N,s] = size(ctrs);
427  DM = zeros(M,N);
428  % Accumulate sum of squares of coordinate differences
429  % The ndgrid command produces two MxN matrices:
430  %   dr, consisting of N identical columns (each containing
431  %       the d-th coordinate of the M data sites)
432  %   cc, consisting of M identical rows (each containing
433  %       the d-th coordinate of the N centers)
434  for d=1:s
435     [dr,cc] = ndgrid(dsites(:,d),ctrs(:,d));
436     DM = DM + (dr-cc).^2;
437  end
438  DM = sqrt(DM);
439
440
441  % rbf = tps(e,r)
442% Defines thin plate spline RBF
443function rbf = tps(e,r)
444rbf = zeros(size(r));
445nz = find(r~=0);   % to deal with singularity at origin
446rbf(nz) = (e*r(nz)).^2.*log(e*r(nz));
447
Note: See TracBrowser for help on using the repository browser.