source: trunk/src/civ_uvmat.m @ 242

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

phys_XYZ introduced as a main function, geometry_calib corrected to deal with reverse x axis in mode linear

File size: 18.7 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    [Ures, Vres,FFres]=...
132                            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);
133                        size(Ures)
134                        size(Vres)
135                        size(FFres)
136                        size(ind_good)
137      Data.Patch1_U(ind_good)=Ures;
138      Data.Patch1_V(ind_good)=Vres;
139      Data.Civ1_FF(ind_good)=FFres
140      Data.CivStage=3;                               
141end   
142%% write result
143errormsg=struct2nc(ncfile,Data);
144
145
146
147%'RUN_FIX': function for fixing velocity fields:
148%-----------------------------------------------
149% RUN_FIX(filename,field,flagindex,thresh_vecC,thresh_vel,iter,flag_mask,maskname,fileref,fieldref)
150%
151%filename: name of the netcdf file (used as input and output)
152%field: structure specifying the names of the fields to fix (depending on civ1 or civ2)
153    %.vel_type='civ1' or 'civ2';
154    %.nb=name of the dimension common to the field to fix ('nb_vectors' for civ1);
155    %.fixflag=name of fix flag variable ('vec_FixFlag' for civ1)
156%flagindex: flag specifying which values of vec_f are removed:
157        % if flagindex(1)=1: vec_f=-2 vectors are removed
158        % if flagindex(2)=1: vec_f=3 vectors are removed
159        % if flagindex(3)=1: vec_f=2 vectors are removed (if iter=1) or vec_f=4 vectors are removed (if iter=2)
160%iter=1 for civ1 fields and iter=2 for civ2 fields
161%thresh_vecC: threshold in the image correlation vec_C
162%flag_mask: =1 mask used to remove vectors (0 else)
163%maskname: name of the mask image file for fix
164%thresh_vel: threshold on velocity, or on the difference with the reference file fileref if exists
165%inf_sup=1: remove values smaller than threshold thresh_vel, =2, larger than threshold
166%fileref: .nc file name for a reference velocity (='': refrence 0 used)
167%fieldref: 'civ1','filter1'...feld used in fileref
168
169function FF=fix_uvmat(Param,F,C,U,V)
170%error=[]; %default
171FF=zeros(size(F));%default
172
173%criterium on warn flags
174if isfield (Param,'WarnFlags')
175    for iflag=1:numel(Param.WarnFlags)
176        FF=(FF==1| F==Param.WarnFlags(iflag));
177    end
178end
179
180%criterium on correlation values
181if isfield (Param,'LowerBoundCorr')
182    FF=FF==1 | C<Param.LowerBoundCorr;
183end
184
185if isfield (Param,'LowerBoundVel')&& ~isequal(Param.LowerBoundVel,0)
186    thresh=Param.LowerBoundVel*Param.LowerBoundVel;
187    FF=FF==1 | (U.*U+V.*V)<thresh;
188end
189if isfield (Param,'UpperBoundVel')&& ~isequal(Param.UpperBoundVel,0)
190    thresh=Param.UpperBoundVel*Param.UpperBoundVel;
191    FF=FF==1 | (U.*U+V.*V)>thresh;
192end
193 
194%
195% % criterium on velocity values
196% delta_u=Field.U;%default without ref file
197% delta_v=Field.V;
198% if exist('fileref','var') && ~isempty(fileref)
199%     if ~exist(fileref,'file')
200%         error='reference file not found in RUN_FIX.m';
201%         display(error);
202%         return
203%     end
204%     FieldRef=read_civxdata(fileref,[],fieldref);   
205%     if isfield(FieldRef,'FF')
206%         index_true=find(FieldRef.FF==0);
207%         FieldRef.X=FieldRef.X(index_true);
208%         FieldRef.Y=FieldRef.Y(index_true);
209%         FieldRef.U=FieldRef.U(index_true);
210%         FieldRef.V=FieldRef.V(index_true);
211%     end
212%     if ~isfield(FieldRef,'X') || ~isfield(FieldRef,'Y') || ~isfield(FieldRef,'U') || ~isfield(FieldRef,'V')
213%         error='reference file is not a velocity field in RUN_FIX.m '; %bad input file
214%         return
215%     end
216%     if length(FieldRef.X)<=1
217%         errordlg('reference field with one vector or less in RUN_FIX.m')
218%         return
219%     end
220%     vec_U_ref=griddata_uvmat(FieldRef.X,FieldRef.Y,FieldRef.U,Field.X,Field.Y);  %interpolate vectors in the ref field
221%     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     
222%     delta_u=Field.U-vec_U_ref;%take the difference with the interpolated ref field
223%     delta_v=Field.V-vec_V_ref;
224% end
225% thresh_vel_x=thresh_vel;
226% thresh_vel_y=thresh_vel;
227% if isequal(inf_sup,1)
228%     flag5=abs(delta_u)<thresh_vel_x & abs(delta_v)<thresh_vel_y &(flag1~=1)&(flag2~=1)&(flag3~=1)&(flag4~=1);
229% elseif isequal(inf_sup,2)
230%     flag5=(abs(delta_u)>thresh_vel_x | abs(delta_v)>thresh_vel_y) &(flag1~=1)&(flag2~=1)&(flag3~=1)&(flag4~=1);
231% end
232%
233%             % flag7 introduce a grey mask, matrix M
234% if isequal (flag_mask,1)
235%    M=imread(maskname);
236%    nxy=size(M);
237%    M=reshape(M,1,nxy(1)*nxy(2));
238%    rangx0=[0.5 nxy(2)-0.5];
239%    rangy0=[0.5 nxy(1)-0.5];
240%    vec_x1=Field.X-Field.U/2;%beginning points
241%    vec_x2=Field.X+Field.U/2;%end points of vectors
242%    vec_y1=Field.Y-Field.V/2;%beginning points
243%    vec_y2=Field.Y+Field.V/2;%end points of vectors
244%    indx=1+round((nxy(2)-1)*(vec_x1-rangx0(1))/(rangx0(2)-rangx0(1)));% image index x at abcissa vec_x
245%    indy=1+round((nxy(1)-1)*(vec_y1-rangy0(1))/(rangy0(2)-rangy0(1)));% image index y at ordinate vec_y   
246%    test_in=~(indx < 1 |indy < 1 | indx > nxy(2) |indy > nxy(1)); %=0 out of the mask image, 1 inside
247%    indx=indx.*test_in+(1-test_in); %replace indx by 1 out of the mask range
248%    indy=indy.*test_in+(1-test_in); %replace indy by 1 out of the mask range
249%    ICOMB=((indx-1)*nxy(1)+(nxy(1)+1-indy));%determine the indices in the image reshaped in a Matlab vector
250%    Mvalues=M(ICOMB);
251%    flag7b=((20 < Mvalues) & (Mvalues < 200))| ~test_in';
252%    indx=1+round((nxy(2)-1)*(vec_x2-rangx0(1))/(rangx0(2)-rangx0(1)));% image index x at abcissa Field.X
253%    indy=1+round((nxy(1)-1)*(vec_y2-rangy0(1))/(rangy0(2)-rangy0(1)));% image index y at ordinate vec_y
254%    test_in=~(indx < 1 |indy < 1 | indx > nxy(2) |indy > nxy(1)); %=0 out of the mask image, 1 inside
255%    indx=indx.*test_in+(1-test_in); %replace indx by 1 out of the mask range
256%    indy=indy.*test_in+(1-test_in); %replace indy by 1 out of the mask range
257%    ICOMB=((indx-1)*nxy(1)+(nxy(1)+1-indy));%determine the indices in the image reshaped in a Matlab vector
258%    Mvalues=M(ICOMB);
259%    flag7e=((Mvalues > 20) & (Mvalues < 200))| ~test_in';
260%    flag7=(flag7b|flag7e)';
261% else
262%    flag7=0;
263% end   
264% flagmagenta=flag1|flag2|flag3|flag4|flag5|flag7;
265% fixflag_unit=Field.FF-10*floor(Field.FF/10); %unity term of fix_flag
266
267
268
269%------------------------------------------------------------------------
270% patch function
271function [U_smooth,V_smooth,FF] =patch_uvmat(X,Y,U,V,Rho,Threshold,SubDomain)
272%subdomain decomposition
273warning off
274nbvec=numel(X);
275NbSubDomain=ceil(nbvec/SubDomain);
276MinX=min(X);
277MinY=min(Y);
278MaxX=max(X);
279MaxY=max(Y);
280RangX=MaxX-MinX;
281RangY=MaxY-MinY;
282AspectRatio=RangY/RangX;
283NbSubDomainX=ceil(sqrt(NbSubDomain/AspectRatio))
284NbSubDomainY=ceil(sqrt(NbSubDomain*AspectRatio))
285
286SizX=RangX/NbSubDomainX;%width of subdomains
287SizY=RangY/NbSubDomainY;%height of subdomains
288CentreX=linspace(MinX+SizX/2,MaxX-SizX/2,NbSubDomainX);
289CentreY=linspace(MinY+SizY/2,MaxY-SizY/2,NbSubDomainY);
290SubIndexX=ceil((X-MinX)/SizX);%subdomain index of vectors
291SubIndexY=ceil((Y-MinY)/SizY);
292rho=RangX*RangY*Rho;%optimum rho increase as teh area of the subdomain
293U_tps=zeros(length(X)+3,1);%default spline
294V_tps=zeros(length(X)+3,1);%default spline
295for isubx=1:NbSubDomainX
296    for isuby=1:NbSubDomainY
297        SubRangx(isuby,isubx,:)=[CentreX(isubx)-SizX/2 CentreX(isubx)+SizX/2]
298        SubRangy(isuby,isubx,:)=[CentreY(isubx)-SizY/2 CentreY(isubx)+SizY/2]
299        for iter=1:3 %increase the subdomain during three iterations at most
300            ind_sel=find(X>SubRangx(isuby,isubx,1) & X<SubRangx(isuby,isubx,2) & Y>SubRangy(isuby,isubx,1) & Y<SubRangy(isuby,isubx,2)); 
301            [U_smooth(ind_sel)]=tps_uvmat(X(ind_sel),Y(ind_sel),U(ind_sel),rho);
302            [V_smooth(ind_sel)]=tps_uvmat(X(ind_sel),Y(ind_sel),V(ind_sel),rho);
303            iter
304            numel(ind_sel)
305            if numel(ind_sel)<SubDomain/4;% too few selected vectors, increase the subrange for next iteration
306                SubRangx(isuby,isubx,1)=SubRangx(isuby,isubx,1)-SizX/4;
307                SubRangx(isuby,isubx,2)=SubRangx(isuby,isubx,2)+SizX/4;
308                SubRangy(isuby,isubx,1)=SubRangy(isuby,isubx,1)-SizY/4;
309                SubRangy(isuby,isubx,2)=SubRangy(isuby,isubx,2)+SizY/4;
310            else
311                UDiff=U_smooth(ind_sel)-U(ind_sel);
312                VDiff=U_smooth(ind_sel)-U(ind_sel);
313                NormDiff=UDiff.*UDiff+VDiff.*VDiff;
314                FF(ind_sel)=20*(NormDiff>Threshold);%put FF value to 20 to identify the criterium of elimmination
315                ind_ind_false=find(FF(ind_sel)~=0);
316                ind_ind_sel=find(FF(ind_sel)==0);
317                U_smooth(ind_sel(ind_ind_false))=zeros(numel(ind_ind_false),1);
318                V_smooth(ind_sel(ind_ind_false))=zeros(numel(ind_ind_false),1);         
319                [U_smooth(ind_sel(ind_ind_sel))]=tps_uvmat(X(ind_sel(ind_ind_sel)),Y(ind_sel(ind_ind_sel)),U(ind_sel(ind_ind_sel)),rho);
320                [V_smooth(ind_sel(ind_ind_sel))]=tps_uvmat(X(ind_sel(ind_ind_sel)),Y(ind_sel(ind_ind_sel)),V(ind_sel(ind_ind_sel)),rho);
321                 if numel(ind_ind_sel)<SubDomain/4;% too few selected vectors, increase the subrange for next iteration
322                    SubRangx(isuby,isubx,1)=SubRangx(isuby,isubx,1)-SizX/4;
323                    SubRangx(isuby,isubx,2)=SubRangx(isuby,isubx,2)+SizX/4;
324                    SubRangy(isuby,isubx,1)=SubRangy(isuby,isubx,1)-SizY/4;
325                    SubRangy(isuby,isubx,2)=SubRangy(isuby,isubx,2)+SizY/4;
326                 else
327                     break
328                 end
329            end       
330        end
331        %intrpolate to new points XI, YI
332%         X_ctrs{isuby,isubx}=X(ind_sel(ind_ind_sel));%positions of tps sources for the subdomain i,j
333%         Y_ctrs{isuby,isubx}=Y(ind_sel(ind_ind_sel));     
334    end
335end
336
337
338%------------------------------------------------------------------------
339%fasshauer@iit.edu MATH 590 ? Chapter 19 32
340% X,Y initial coordiantes
341% XI vector, YI column vector for the grid of interpolation points
342function [U_smooth,U_tps]=tps_uvmat(X,Y,U,rho)
343%------------------------------------------------------------------------
344%rho smoothing parameter
345ep = 1;
346X=reshape(X,[],1);
347Y=reshape(Y,[],1);
348rhs = reshape(U,[],1);
349% if exist('FF','var')
350% test_false=isnan(rhs)|FF~=0;
351% else
352%     test_false=isnan(rhs);
353% end
354% X(test_false)=[];
355% Y(test_false)=[];
356% rhs(test_false)=[];
357%randn('state',3); rhs = rhs + 0.03*randn(size(rhs));
358rhs = [rhs; zeros(3,1)];
359dsites = [X Y];% coordinates of measurement sites
360ctrs = dsites;%radial base functions are located at the measurement sites
361DM_data = DistanceMatrix(dsites,ctrs);%2D matrix of distances between spline centres (=initial points) ctrs
362% if size(XI,1)==1 && size(YI,2)==1 % XI vector, YI column vector
363%      [XI,YI]=meshgrid(XI,YI);
364% end
365% [npy,npx]=size(XI);
366% epoints = [reshape(XI,[],1) reshape(YI,[],1)];
367IM_sites = tps(ep,DM_data);%values of thin plate at site points
368IM = IM_sites + rho*eye(size(IM_sites));%  rho=1/(2*omega) , omega given by fasshauer;
369PM=[ones(size(dsites,1),1) dsites];
370IM=[IM PM; [PM' zeros(3,3)]];
371%fprintf('Condition number estimate: %e\n',condest(IM))
372%DM_eval = DistanceMatrix(epoints,ctrs);%2D matrix of distances between extrapolation points epoints and spline centres (=site points) ctrs
373%EM = tps(ep,DM_eval);%values of thin plate
374%PM = [ones(size(epoints,1),1) epoints];
375%EM = [EM PM];
376U_tps=(IM\rhs);
377PM = [ones(size(dsites,1),1) dsites];
378EM = [IM_sites PM];
379U_smooth=EM *U_tps;
380
381
382% U_patch = EM * spline_coeff;
383% U_patch=reshape(U_patch,npy,npx);
384% PM = [ones(size(dsites,1),1) dsites];
385% EM = [IM_sites PM];
386% U(test_false)=[];
387% U_nodes=EM * spline_coeff;
388
389%exact = testfunctions(epoints);
390%maxerr = norm(Pf-exact,inf);
391% PlotSurf(xe,ye,Pf,neval,exact,maxerr,[160,20]);
392% PlotError2D(xe,ye,Pf,exact,maxerr,neval,[160,20]);
393
394
395
396  % DM = DistanceMatrix(dsites,ctrs)
397% Forms the distance matrix of two sets of points in R^s,
398% i.e., DM(i,j) = || datasite_i - center_j ||_2.
399% Input
400%   dsites: Mxs matrix representing a set of M data sites in R^s
401%              (i.e., each row contains one s-dimensional point)
402%   ctrs:   Nxs matrix representing a set of N centers in R^s
403%              (one center per row)
404% Output
405
406
407%   DM:     MxN matrix whose i,j position contains the Euclidean
408%              distance between the i-th data site and j-th center
409  function DM = DistanceMatrix(dsites,ctrs)
410  [M,s] = size(dsites); [N,s] = size(ctrs);
411  DM = zeros(M,N);
412  % Accumulate sum of squares of coordinate differences
413  % The ndgrid command produces two MxN matrices:
414  %   dr, consisting of N identical columns (each containing
415  %       the d-th coordinate of the M data sites)
416  %   cc, consisting of M identical rows (each containing
417  %       the d-th coordinate of the N centers)
418  for d=1:s
419     [dr,cc] = ndgrid(dsites(:,d),ctrs(:,d));
420     DM = DM + (dr-cc).^2;
421  end
422  DM = sqrt(DM);
423
424
425  % rbf = tps(e,r)
426% Defines thin plate spline RBF
427function rbf = tps(e,r)
428rbf = zeros(size(r));
429nz = find(r~=0);   % to deal with singularity at origin
430rbf(nz) = (e*r(nz)).^2.*log(e*r(nz));
431
Note: See TracBrowser for help on using the repository browser.