source: trunk/src/civ_matlab.m @ 388

Last change on this file since 388 was 388, checked in by sommeria, 12 years ago

several bugs corrected
file indexing, color images...

File size: 31.2 KB
Line 
1%'civ_matlab': Matlab version of the PIV programs CivX
2% --- call the sub-functions:
3%   civ: PIV function itself
4%   fix: removes false vectors after detection by various criteria
5%   filter_tps: make interpolation-smoothing
6%------------------------------------------------------------------------
7% function [Data,errormsg,result_conv]= civ_uvmat(Param,ncfile)
8%
9%OUTPUT
10% Data=structure containing the PIV results and information on the processing parameters
11% errormsg=error message char string, default=''
12% resul_conv: image inter-correlation function for the last grid point (used for tests)
13%
14%INPUT:
15% Param: input images and processing parameters
16%     .Civ1: for civ1
17%     .Fix1:
18%     .Patch1:
19%     .Civ2: for civ2
20%     .Fix2:
21%     .Patch2:
22% ncfile: name of a netcdf file to be created for the result (extension .nc)
23%
24%AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA
25%  Copyright  2011, LEGI / CNRS-UJF-INPG, joel.sommeria@legi.grenoble-inp.fr.
26%AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA
27%     This is part of the toolbox UVMAT.
28%
29%     UVMAT is free software; you can redistribute it and/or modify
30%     it under the terms of the GNU General Public License as published by
31%     the Free Software Foundation; either version 2 of the License, or
32%     (at your option) any later version.
33%
34%     UVMAT is distributed in the hope that it will be useful,
35%     but WITHOUT ANY WARRANTY; without even the implied warranty of
36%     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
37%     GNU General Public License (open UVMAT/COPYING.txt) for more details.
38%AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA
39
40function [Data,errormsg,result_conv]= civ_matlab(Param,ncfile)
41errormsg='';
42Data.ListGlobalAttribute={'Conventions','Program','CivStage'};
43Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
44Data.Program='civ_matlab';
45Data.CivStage=0;%default
46ListVarCiv1={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F'}; %variables to read
47ListVarFix1={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F','Civ1_FF'};
48mask='';
49maskname='';%default
50check_civx=0;%default
51check_civ1=0;%default
52check_patch1=0;%default
53
54if ischar(Param)
55    Param=xml2struct(Param); %if Param is the name of an xml file, read this file as a Matlab structure
56end
57
58%% Civ1
59if isfield (Param,'Civ1')
60    check_civ1=1;% test for further use of civ1 results
61    % %% prepare images
62    par_civ1=Param.Civ1;
63    if isfield(par_civ1,'reverse_pair')
64        if par_civ1.reverse_pair
65            if ischar(par_civ1.ImageB)
66                temp=par_civ1.ImageA;
67                par_civ1.ImageA=imread(par_civ1.ImageB);
68            end
69            if ischar(temp)
70                par_civ1.ImageB=imread(temp);
71            end
72        end
73    else
74        if isfield(par_civ1,'ImageA') && ischar(par_civ1.ImageA) % case with no image: only the PIV grid is calculated
75            [Field,ParamOut,errormsg] = read_field(par_civ1.ImageA,par_civ1.FileTypeA,[],par_civ1.i1);
76            par_civ1.ImageA=Field.A;%imread(par_civ1.ImageA);%[Field,ParamOut,errormsg] = read_field(ObjectName,FileType,ParamIn,num)
77        end
78        if isfield(par_civ1,'ImageB')&& ischar(par_civ1.ImageB)
79            [Field,ParamOut,errormsg] = read_field(par_civ1.ImageB,par_civ1.FileTypeB,[],par_civ1.i2);
80            par_civ1.ImageB=Field.A;%=imread(par_civ1.ImageB);
81        end
82    end
83   
84    % caluclate velocity data (y and v in indices, reverse to y component)
85    [xtable ytable utable vtable ctable F result_conv errormsg] = civ (par_civ1);
86   
87    % to try the reverse_pair method, uncomment below
88    %     [xtable1 ytable1 utable1 vtable1 ctable1 F1 result_conv1 errormsg1] = civ (Param.Civ1);
89    %     Param.Civ1.reverse_pair=1;
90    %     [xtable2 ytable2 utable2 vtable2 ctable2 F2 result_conv2 errormsg2] = civ (Param.Civ1);
91    %     xtable=[xtable1; xtable2];
92    %     ytable=[ytable1; ytable2];
93    %     utable=[utable1; -utable2];
94    %     vtable=[vtable1; -vtable2];
95    %     ctable=[ctable1; ctable2];
96    %     F=[F1; F2];
97    %     result_conv=[result_conv1; result_conv2];
98    %     errormsg=[errormsg1; errormsg2];
99    if ~isempty(errormsg)
100        return
101    end
102    list_param=(fieldnames(Param.Civ1))';
103    Civ1_param=list_param;%default
104    for ilist=1:length(list_param)
105        Civ1_param{ilist}=['Civ1_' list_param{ilist}];
106        Data.(['Civ1_' list_param{ilist}])=Param.Civ1.(list_param{ilist});
107    end
108    Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ1_param];% {'Civ1_Time','Civ1_Dt'}];
109    Data.ListVarName={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_F','Civ1_C'};%  cell array containing the names of the fields to record
110    Data.VarDimName={'NbVec1','NbVec1','NbVec1','NbVec1','NbVec1','NbVec1'};
111    Data.VarAttribute{1}.Role='coord_x';
112    Data.VarAttribute{2}.Role='coord_y';
113    Data.VarAttribute{3}.Role='vector_x';
114    Data.VarAttribute{4}.Role='vector_y';
115    Data.VarAttribute{5}.Role='warnflag';
116    Data.Civ1_X=reshape(xtable,[],1);
117    Data.Civ1_Y=reshape(Param.Civ1.ImageHeight-ytable+1,[],1);
118    Data.Civ1_U=reshape(utable,[],1);
119    Data.Civ1_V=reshape(-vtable,[],1);
120    Data.Civ1_C=reshape(ctable,[],1);
121    Data.Civ1_F=reshape(F,[],1);
122    Data.CivStage=1; 
123else
124    if exist('ncfile','var')
125        CivFile=ncfile;
126    elseif isfield(Param.Patch1,'CivFile')
127        CivFile=Param.Patch1.CivFile;
128    end
129    Data=nc2struct(CivFile,'ListGlobalAttribute','absolut_time_T0'); %look for the constant 'absolut_time_T0' to detect old civx data format
130    if isfield(Data,'Txt')
131        errormsg=Data.Txt;
132        return
133    end
134    if ~isempty(Data.absolut_time_T0')%read civx file
135        check_civx=1;% test for old civx data format
136        [Data,vardetect,ichoice]=nc2struct(CivFile);%read the variables in the netcdf file
137    else
138        Data=nc2struct(CivFile);%read civ1 and fix1 data in the existing netcdf file
139    end
140end
141
142%% Fix1
143if isfield (Param,'Fix1')
144    ListFixParam=fieldnames(Param.Fix1);
145    for ilist=1:length(ListFixParam)
146        ParamName=ListFixParam{ilist};
147        ListName=['Fix1_' ParamName];
148        eval(['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];'])
149        eval(['Data.' ListName '=Param.Fix1.' ParamName ';'])
150    end
151    if check_civx
152        if ~isfield(Data,'fix')
153            Data.ListGlobalAttribute=[Data.ListGlobalAttribute 'fix'];
154            Data.fix=1;
155            Data.ListVarName=[Data.ListVarName {'vec_FixFlag'}];
156            Data.VarDimName=[Data.VarDimName {'nb_vectors'}];
157        end
158        Data.vec_FixFlag=fix(Param.Fix1,Data.vec_F,Data.vec_C,Data.vec_U,Data.vec_V,Data.vec_X,Data.vec_Y);
159    else
160        Data.ListVarName=[Data.ListVarName {'Civ1_FF'}];
161        Data.VarDimName=[Data.VarDimName {'NbVec1'}];
162        nbvar=length(Data.ListVarName);
163        Data.VarAttribute{nbvar}.Role='errorflag';   
164        Data.Civ1_FF=fix(Param.Fix1,Data.Civ1_F,Data.Civ1_C,Data.Civ1_U,Data.Civ1_V);
165        Data.CivStage=2;   
166    end
167end   
168%% Patch1
169if isfield (Param,'Patch1')
170    if check_civx
171        errormsg='Civ Matlab input needed for patch';
172        return
173    end
174    check_patch1=1;
175    Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Patch1_Rho','Patch1_Threshold','Patch1_SubDomain'}];
176    Data.Patch1_Rho=Param.Patch1.SmoothingParam;
177    Data.Patch1_Threshold=Param.Patch1.MaxDiff;
178    Data.Patch1_SubDomain=Param.Patch1.SubdomainSize;
179    Data.ListVarName=[Data.ListVarName {'Civ1_U_Diff','Civ1_V_Diff','Civ1_SubRange','Civ1_NbSites','Civ1_Coord_tps','Civ1_U_tps','Civ1_V_tps'}];
180    Data.VarDimName=[Data.VarDimName {'NbVec1','NbVec1',{'NbCoord','Two','NbSubDomain1'},{'NbSubDomain1'},...
181             {'NbVec1Sub','NbCoord','NbSubDomain1'},{'Nbtps1','NbSubDomain1'},{'Nbtps1','NbSubDomain1'}}];
182    nbvar=length(Data.ListVarName);
183    Data.VarAttribute{nbvar-1}.Role='vector_x';
184    Data.VarAttribute{nbvar}.Role='vector_y';
185    Data.Civ1_U_Diff=zeros(size(Data.Civ1_X));
186    Data.Civ1_V_Diff=zeros(size(Data.Civ1_X));
187    if isfield(Data,'Civ1_FF')
188        ind_good=find(Data.Civ1_FF==0);
189    else
190        ind_good=1:numel(Data.Civ1_X);
191    end
192    [Data.Civ1_SubRange,Data.Civ1_NbSites,Data.Civ1_Coord_tps,Data.Civ1_U_tps,Data.Civ1_V_tps,tild,Ures, Vres,tild,FFres]=...
193            filter_tps([Data.Civ1_X(ind_good) Data.Civ1_Y(ind_good)],Data.Civ1_U(ind_good),Data.Civ1_V(ind_good),[],Data.Patch1_SubDomain,Data.Patch1_Rho,Data.Patch1_Threshold);
194      Data.Civ1_U_Diff(ind_good)=Data.Civ1_U(ind_good)-Ures;
195      Data.Civ1_V_Diff(ind_good)=Data.Civ1_V(ind_good)-Vres;
196      Data.Civ1_FF(ind_good)=FFres;
197      Data.CivStage=3;                             
198end   
199
200%% Civ2
201if isfield (Param,'Civ2')
202    par_civ2=Param.Civ2;
203    if ~isfield (Param,'Civ1') || ~strcmp(Param.Civ1.ImageA,par_civ2.ImageA)
204        par_civ2.ImageA=imread(Param.Civ2.ImageA);%read first image if not already done for civ1
205    else
206        par_civ2.ImageA=par_civ1.ImageA;
207    end
208    if ~isfield (Param,'Civ1') || ~strcmp(Param.Civ1.ImageB,par_civ2.ImageB)
209        par_civ2.ImageB=imread(Param.Civ2.ImageB);%read second image if not already done for civ1
210         else
211        par_civ2.ImageB=par_civ1.ImageB;
212    end
213    ibx2=ceil(par_civ2.Bx/2);
214    iby2=ceil(par_civ2.By/2);
215    isx2=ibx2+5;% search ara +-5 pixels around the guess
216    isy2=iby2+5;
217    % shift from par_civ2.filename_nc1
218    % shiftx=velocity interpolated at position
219    miniy=max(1+isy2,1+iby2);
220    minix=max(1+isx2,1+ibx2);
221    maxiy=min(size(par_civ2.ImageA,1)-isy2,size(par_civ2.ImageA,1)-iby2);
222    maxix=min(size(par_civ2.ImageA,2)-isx2,size(par_civ2.ImageA,2)-ibx2);
223    [GridX,GridY]=meshgrid(minix:par_civ2.Dx:maxix,miniy:par_civ2.Dy:maxiy);
224    GridX=reshape(GridX,[],1);
225    GridY=reshape(GridY,[],1);
226    Shiftx=zeros(size(GridX));% shift expected from civ1 data
227    Shifty=zeros(size(GridX));
228    nbval=zeros(size(GridX));
229    if par_civ2.CheckDeformation
230        DUDX=zeros(size(GridX));
231        DUDY=zeros(size(GridX));
232        DVDX=zeros(size(GridX));
233        DVDY=zeros(size(GridX));
234    end
235    NbSubDomain=size(Data.Civ1_SubRange,3);
236    % get the guess from patch1
237    for isub=1:NbSubDomain
238        nbvec_sub=Data.Civ1_NbSites(isub);
239        ind_sel=find(GridX>=Data.Civ1_SubRange(1,1,isub) & GridX<=Data.Civ1_SubRange(1,2,isub) & GridY>=Data.Civ1_SubRange(2,1,isub) & GridY<=Data.Civ1_SubRange(2,2,isub));
240        epoints = [GridX(ind_sel) GridY(ind_sel)];% coordinates of interpolation sites
241        ctrs=Data.Civ1_Coord_tps(1:nbvec_sub,:,isub) ;%(=initial points) ctrs
242        nbval(ind_sel)=nbval(ind_sel)+1;% records the number of values for eacn interpolation point (in case of subdomain overlap)
243        EM = tps_eval(epoints,ctrs);
244        Shiftx(ind_sel)=Shiftx(ind_sel)+EM*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
245        Shifty(ind_sel)=Shifty(ind_sel)+EM*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
246        if par_civ2.CheckDeformation
247            [EMDX,EMDY] = tps_eval_dxy(epoints,ctrs);%2D matrix of distances between extrapolation points epoints and spline centres (=site points) ctrs
248            DUDX(ind_sel)=DUDX(ind_sel)+EMDX*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
249            DUDY(ind_sel)=DUDY(ind_sel)+EMDY*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
250            DVDX(ind_sel)=DVDX(ind_sel)+EMDX*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
251            DVDY(ind_sel)=DVDY(ind_sel)+EMDY*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
252        end
253    end
254    mask='';
255    if par_civ2.CheckMask&&~isempty(par_civ2.Mask)&& ~strcmp(maskname,par_civ2.Mask)% mask exist, not already read in civ1
256        mask=imread(par_civ2.Mask);
257    end
258    par_civ2.Searchx=2*isx2+1;
259    par_civ2.Searchy=2*isy2+1;
260    par_civ2.Shiftx=Shiftx(nbval>=1)./nbval(nbval>=1);
261    par_civ2.Shifty=Shifty(nbval>=1)./nbval(nbval>=1);
262    par_civ2.Grid=[GridX(nbval>=1)-par_civ2.Shiftx/2 GridY(nbval>=1)-par_civ2.Shifty/2];% grid taken at the extrapolated origin of the displacement vectors   
263    if par_civ2.CheckDeformation
264        par_civ2.DUDX=DUDX./nbval;
265        par_civ2.DUDY=DUDY./nbval;
266        par_civ2.DVDX=DVDX./nbval;
267        par_civ2.DVDY=DVDY./nbval;
268    end
269    % caluclate velocity data (y and v in indices, reverse to y component)
270    [xtable ytable utable vtable ctable F] = civ (par_civ2);
271%     diff_squared=(utable-par_civ2.Shiftx).*(utable-par_civ2.Shiftx)+(vtable+par_civ2.Shifty).*(vtable+par_civ2.Shifty);
272%     F(diff_squared>=4)=4; %flag vectors whose distance to the guess exceeds 2 pixels
273    list_param=(fieldnames(Param.Civ2))';
274    list_remove={'pxcmx','pxcmy','npx','npy','gridflag','maskflag','term_a','term_b','T0'};
275    for ilist=1:length(list_remove)
276        index=strcmp(list_remove{ilist},list_param);
277        if ~isempty(find(index,1))
278            list_param(index)=[];
279        end
280    end
281    for ilist=1:length(list_param)
282        Civ2_param{ilist}=['Civ2_' list_param{ilist}];
283        eval(['Data.Civ2_' list_param{ilist} '=Param.Civ2.' list_param{ilist} ';'])
284    end
285    if isfield(Data,'Civ2_gridname') && strcmp(Data.Civ1_gridname(1:6),'noFile')
286        Data.Civ1_gridname='';
287    end
288    if isfield(Data,'Civ2_maskname') && strcmp(Data.Civ1_maskname(1:6),'noFile')
289        Data.Civ2_maskname='';
290    end
291    Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ2_param {'Civ2_Time','Civ2_Dt'}];
292%     Data.Civ2_Time=par_civ2.Time;
293%     Data.Civ2_Dt=par_civ2.Dt;
294    nbvar=numel(Data.ListVarName);
295    Data.ListVarName=[Data.ListVarName {'Civ2_X','Civ2_Y','Civ2_U','Civ2_V','Civ2_F','Civ2_C'}];%  cell array containing the names of the fields to record
296    Data.VarDimName=[Data.VarDimName {'NbVec2','NbVec2','NbVec2','NbVec2','NbVec2','NbVec2'}];
297    Data.VarAttribute{nbvar+1}.Role='coord_x';
298    Data.VarAttribute{nbvar+2}.Role='coord_y';
299    Data.VarAttribute{nbvar+3}.Role='vector_x';
300    Data.VarAttribute{nbvar+4}.Role='vector_y';
301    Data.VarAttribute{nbvar+5}.Role='warnflag';
302    Data.Civ2_X=reshape(xtable,[],1);
303    Data.Civ2_Y=reshape(size(par_civ2.ImageA,1)-ytable+1,[],1);
304    Data.Civ2_U=reshape(utable,[],1);
305    Data.Civ2_V=reshape(-vtable,[],1);
306    Data.Civ2_C=reshape(ctable,[],1);
307    Data.Civ2_F=reshape(F,[],1);
308    Data.CivStage=Data.CivStage+1;
309end
310
311%% Fix2
312if isfield (Param,'Fix2')
313    ListFixParam=fieldnames(Param.Fix2);
314    for ilist=1:length(ListFixParam)
315        ParamName=ListFixParam{ilist};
316        ListName=['Fix2_' ParamName];
317        eval(['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];'])
318        eval(['Data.' ListName '=Param.Fix2.' ParamName ';'])
319    end
320    if check_civx
321        if ~isfield(Data,'fix2')
322            Data.ListGlobalAttribute=[Data.ListGlobalAttribute 'fix2'];
323            Data.fix2=1;
324            Data.ListVarName=[Data.ListVarName {'vec2_FixFlag'}];
325            Data.VarDimName=[Data.VarDimName {'nb_vectors2'}];
326        end
327        Data.vec_FixFlag=fix(Param.Fix2,Data.vec2_F,Data.vec2_C,Data.vec2_U,Data.vec2_V,Data.vec2_X,Data.vec2_Y);
328    else
329        Data.ListVarName=[Data.ListVarName {'Civ2_FF'}];
330        Data.VarDimName=[Data.VarDimName {'nbvec2'}];
331        nbvar=length(Data.ListVarName);
332        Data.VarAttribute{nbvar}.Role='errorflag';   
333        Data.Civ2_FF=fix(Param.Fix2,Data.Civ2_F,Data.Civ2_C,Data.Civ2_U,Data.Civ2_V);
334        Data.CivStage=Data.CivStage+1;   
335    end
336   
337end   
338
339%% Patch2
340if isfield (Param,'Patch2')
341    Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Patch2_Rho','Patch2_Threshold','Patch2_SubDomain'}];
342    Data.Patch2_Rho=Param.Patch2.SmoothingParam;
343    Data.Patch2_Threshold=Param.Patch2.MaxDiff;
344    Data.Patch2_SubDomain=Param.Patch2.SubdomainSize;
345    Data.ListVarName=[Data.ListVarName {'Civ2_U_Diff','Civ2_V_Diff','Civ2_SubRange','Civ2_NbSites','Civ2_Coord_tps','Civ2_U_tps','Civ2_V_tps'}];
346    Data.VarDimName=[Data.VarDimName {'NbVec2','NbVec2',{'NbCoord','Two','NbSubDomain2'},{'NbSubDomain2'},...
347             {'NbVec2Sub','NbCoord','NbSubDomain2'},{'Nbtps2','NbSubDomain2'},{'Nbtps2','NbSubDomain2'}}];
348%     Data.ListVarName=[Data.ListVarName {'Civ2_U_Diff','Civ2_V_Diff','Civ2_X_SubRange','Civ2_Y_SubRange','Civ2_X_tps','Civ2_Y_tps','Civ2_U_tps','Civ2_V_tps','Civ2_Indices_tps'}];
349%     Data.VarDimName=[Data.VarDimName {'NbVec2','NbVec2',{'NbSubDomain2','Two'},{'NbSubDomain2','Two'},...
350%              {'NbVec2Sub','NbSubDomain2'},{'NbVec2Sub','NbSubDomain2'},{'Nbtps2','NbSubDomain2'},{'Nbtps2','NbSubDomain2'},{'NbVec2Sub','NbSubDomain2'}}];
351    nbvar=length(Data.ListVarName);
352    Data.VarAttribute{nbvar-1}.Role='vector_x';
353    Data.VarAttribute{nbvar}.Role='vector_y';
354    Data.Civ2_U_Diff=zeros(size(Data.Civ2_X));
355    Data.Civ2_V_Diff=zeros(size(Data.Civ2_X));
356    if isfield(Data,'Civ2_FF')
357        ind_good=find(Data.Civ2_FF==0);
358    else
359        ind_good=1:numel(Data.Civ2_X);
360    end
361    [Data.Civ2_SubRange,Data.Civ2_NbSites,Data.Civ2_Coord_tps,Data.Civ2_U_tps,Data.Civ2_V_tps,tild,Ures, Vres,tild,FFres]=...
362         filter_tps([Data.Civ2_X(ind_good) Data.Civ2_Y(ind_good)],Data.Civ2_U(ind_good),Data.Civ2_V(ind_good),[],Data.Patch2_SubDomain,Data.Patch2_Rho,Data.Patch2_Threshold);
363      Data.Civ2_U_Diff(ind_good)=Data.Civ2_U(ind_good)-Ures;
364      Data.Civ2_V_Diff(ind_good)=Data.Civ2_V(ind_good)-Vres;
365      Data.Civ2_FF(ind_good)=FFres;
366      Data.CivStage=Data.CivStage+1;                             
367end 
368
369%% write result in a netcdf file if requested
370if exist('ncfile','var')
371    errormsg=struct2nc(ncfile,Data);
372end
373
374% 'civ': function piv.m adapted from PIVlab http://pivlab.blogspot.com/
375%--------------------------------------------------------------------------
376% function [xtable ytable utable vtable typevector] = civ (image1,image2,ibx,iby step, subpixfinder, mask, roi)
377%
378% OUTPUT:
379% xtable: set of x coordinates
380% ytable: set of y coordiantes
381% utable: set of u displacements (along x)
382% vtable: set of v displacements (along y)
383% ctable: max image correlation for each vector
384% typevector: set of flags, =1 for good, =0 for NaN vectors
385%
386%INPUT:
387% image1:first image (matrix)
388% image2: second image (matrix)
389% ibx2,iby2: half size of the correlation box along x and y, in px (size=(2*iby2+1,2*ibx2+1)
390% isx2,isy2: half size of the search box along x and y, in px (size=(2*isy2+1,2*isx2+1)
391% shiftx, shifty: shift of the search box (in pixel index, yshift reversed)
392% step: mesh of the measurement points (in px)
393% subpixfinder=1 or 2 controls the curve fitting of the image correlation
394% mask: =[] for no mask
395% roi: 4 element vector defining a region of interest: x position, y position, width, height, (in image indices), for the whole image, roi=[];
396function [xtable ytable utable vtable ctable F result_conv errormsg] = civ (par_civ)
397%this funtion performs the DCC PIV analysis. Recent window-deformation
398%methods perform better and will maybe be implemented in the future.
399
400%% prepare grid
401ibx2=ceil(par_civ.Bx/2);
402iby2=ceil(par_civ.By/2);
403isx2=ceil(par_civ.Searchx/2);
404isy2=ceil(par_civ.Searchy/2);
405shiftx=round(par_civ.Shiftx);
406shifty=-round(par_civ.Shifty);% sign minus because image j index increases when y decreases
407if isfield(par_civ,'Grid')
408    if ischar(par_civ.Grid)%read the drid file if the input is a file name
409        par_civ.Grid=dlmread(par_civ.Grid);
410        par_civ.Grid(1,:)=[];%the first line must be removed (heading in the grid file)
411    end
412else% automatic measurement grid
413    ibx2=ceil(par_civ.Bx/2);
414    iby2=ceil(par_civ.By/2);
415    isx2=ceil(par_civ.Searchx/2);
416    isy2=ceil(par_civ.Searchy/2);
417    miniy=max(1+isy2+shifty,1+iby2);
418    minix=max(1+isx2-shiftx,1+ibx2);
419    maxiy=min(par_civ.ImageHeight-isy2+shifty,par_civ.ImageHeight-iby2);
420    maxix=min(par_civ.ImageWidth-isx2-shiftx,par_civ.ImageWidth-ibx2);
421    [GridX,GridY]=meshgrid(minix:par_civ.Dx:maxix,miniy:par_civ.Dy:maxiy);
422    par_civ.Grid(:,1)=reshape(GridX,[],1);
423    par_civ.Grid(:,2)=reshape(GridY,[],1);
424end
425nbvec=size(par_civ.Grid,1);
426if numel(shiftx)==1
427    shiftx=shiftx*ones(nbvec,1);
428    shifty=shifty*ones(nbvec,1);
429end
430%% Default output
431xtable=par_civ.Grid(:,1);
432ytable=par_civ.Grid(:,2);
433utable=zeros(nbvec,1);
434vtable=zeros(nbvec,1);
435ctable=zeros(nbvec,1);
436F=zeros(nbvec,1);
437result_conv=[];
438errormsg='';
439
440%% prepare mask
441if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
442    if strcmp(par_civ.Mask,'all')
443        return    % get the grid only, no civ calculation
444    elseif ischar(par_civ.Mask)
445        par_civ.Mask=imread(par_civ.Mask);
446    end
447end
448check_MinIma=isfield(par_civ,'MinIma');% test for image luminosity threshold
449check_MaxIma=isfield(par_civ,'MaxIma') && ~isempty(par_civ.MaxIma);
450
451% %% prepare images
452% if isfield(par_civ,'reverse_pair')
453%     if par_civ.reverse_pair
454%         if ischar(par_civ.ImageB)
455%             temp=par_civ.ImageA;
456%             par_civ.ImageA=imread(par_civ.ImageB);
457%         end
458%         if ischar(temp)
459%             par_civ.ImageB=imread(temp);
460%         end
461%     end
462% else
463%     if ischar(par_civ.ImageA)
464%         par_civ.ImageA=imread(par_civ.ImageA);
465%     end
466%     if ischar(par_civ.ImageB)
467%         par_civ.ImageB=imread(par_civ.ImageB);
468%     end
469% end
470par_civ.ImageA=sum(double(par_civ.ImageA),3);%sum over rgb component for color images
471par_civ.ImageB=sum(double(par_civ.ImageB),3);
472[npy_ima npx_ima]=size(par_civ.ImageA);
473if ~isequal(size(par_civ.ImageB),[npy_ima npx_ima])
474    errormsg='image pair with unequal size';
475    return
476end
477
478%% Apply mask
479    % Convention for mask
480    % mask >200 : velocity calculated
481    %  200 >=mask>150;velocity not calculated, interpolation allowed (bad spots)
482    % 150>=mask >100: velocity not calculated, nor interpolated
483    %  100>=mask> 20: velocity not calculated, impermeable (no flux through mask boundaries) TO IMPLEMENT
484    %  20>=mask: velocity=0
485checkmask=0;
486if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
487   checkmask=1;
488   if ~isequal(size(par_civ.Mask),[npy_ima npx_ima])
489        errormsg='mask must be an image with the same size as the images';
490        return
491   end
492  %  check_noflux=(par_civ.Mask<100) ;%TODO: to implement
493    check_undefined=(par_civ.Mask<200 & par_civ.Mask>=100 );
494    par_civ.ImageA(check_undefined)=min(min(par_civ.ImageA));% put image A to zero (i.e. the min image value) in the undefined  area
495    par_civ.ImageB(check_undefined)=min(min(par_civ.ImageB));% put image B to zero (i.e. the min image value) in the undefined  area
496end
497
498%% compute image correlations: MAINLOOP on velocity vectors
499corrmax=0;
500sum_square=1;% default
501mesh=1;% default
502CheckDecimal=isfield(par_civ,'CheckDecimal')&& par_civ.CheckDecimal==1;
503if CheckDecimal
504    mesh=0.2;%mesh in pixels for subpixel image interpolation
505    CheckDeformation=isfield(par_civ,'CheckDeformation')&& par_civ.CheckDeformation==1;
506end
507% vector=[0 0];%default
508for ivec=1:nbvec
509    iref=round(par_civ.Grid(ivec,1)+0.5);% xindex on the image A for the middle of the correlation box
510    jref=round(par_civ.ImageHeight-par_civ.Grid(ivec,2)+0.5);% yindex on the image B for the middle of the correlation box
511    if ~(checkmask && par_civ.Mask(jref,iref)<=20) %velocity not set to zero by the black mask
512        if jref-iby2<1 || jref+iby2>par_civ.ImageHeight|| iref-ibx2<1 || iref+ibx2>par_civ.ImageWidth||...
513              jref+shifty(ivec)-isy2<1||jref+shifty(ivec)+isy2>par_civ.ImageHeight|| iref+shiftx(ivec)-isx2<1 || iref+shiftx(ivec)+isx2>par_civ.ImageWidth  % we are outside the image
514            F(ivec)=3;
515        else
516            image1_crop=par_civ.ImageA(jref-iby2:jref+iby2,iref-ibx2:iref+ibx2);%extract a subimage (correlation box) from image A
517            image2_crop=par_civ.ImageB(jref+shifty(ivec)-isy2:jref+shifty(ivec)+isy2,iref+shiftx(ivec)-isx2:iref+shiftx(ivec)+isx2);%extract a larger subimage (search box) from image B
518            image1_mean=mean(mean(image1_crop));
519            image2_mean=mean(mean(image2_crop));
520            %threshold on image minimum
521            if check_MinIma && (image1_mean < par_civ.MinIma || image2_mean < par_civ.MinIma)
522                F(ivec)=3;
523            end
524            %threshold on image maximum
525            if check_MaxIma && (image1_mean > par_civ.MaxIma || image2_mean > par_civ.MaxIma)
526                F(ivec)=3;
527            end
528        end     
529        if F(ivec)~=3
530            image1_crop=image1_crop-image1_mean;%substract the mean
531            image2_crop=image2_crop-image2_mean;
532            if CheckDecimal
533                xi=(1:mesh:size(image1_crop,2));
534                yi=(1:mesh:size(image1_crop,1))';
535                if CheckDeformation
536                    [XI,YI]=meshgrid(xi-ceil(size(image1_crop,2)/2),yi-ceil(size(image1_crop,1)/2));
537                    XIant=XI-par_civ.DUDX(ivec)*XI-par_civ.DUDY(ivec)*YI+ceil(size(image1_crop,2)/2);
538                    YIant=YI-par_civ.DVDX(ivec)*XI-par_civ.DVDY(ivec)*YI+ceil(size(image1_crop,1)/2);
539                    image1_crop=interp2(image1_crop,XIant,YIant);
540                else
541                    image1_crop=interp2(image1_crop,xi,yi);
542                end
543                xi=(1:mesh:size(image2_crop,2));
544                yi=(1:mesh:size(image2_crop,1))';
545                image2_crop=interp2(image2_crop,xi,yi);
546            end
547            sum_square=sum(sum(image1_crop.*image1_crop));
548            %reference: Oliver Pust, PIV: Direct Cross-Correlation
549            result_conv= conv2(image2_crop,flipdim(flipdim(image1_crop,2),1),'valid');
550            corrmax= max(max(result_conv));
551            result_conv=(result_conv/corrmax)*255; %normalize, peak=always 255
552            %Find the correlation max, at 255
553            [y,x] = find(result_conv==255,1);
554            if ~isempty(y) && ~isempty(x)
555                try
556                    if par_civ.Rho==1
557                        [vector,F(ivec)] = SUBPIXGAUSS (result_conv,x,y);
558                    elseif par_civ.Rho==2
559                        [vector,F(ivec)] = SUBPIX2DGAUSS (result_conv,x,y);
560                    end
561                    utable(ivec)=vector(1)*mesh+shiftx(ivec);
562                    vtable(ivec)=vector(2)*mesh+shifty(ivec);                 
563                    xtable(ivec)=iref+utable(ivec)/2-0.5;% convec flow (velocity taken at the point middle from imgae 1 and 2)
564                    ytable(ivec)=jref+vtable(ivec)/2-0.5;% and position of pixel 1=0.5 (convention for image coordinates=0 at the edge)
565                    iref=round(xtable(ivec));% image index for the middle of the vector
566                    jref=round(ytable(ivec));
567                    if checkmask && par_civ.Mask(jref,iref)<200 && par_civ.Mask(jref,iref)>=100
568                        utable(ivec)=0;
569                        vtable(ivec)=0;
570                        F(ivec)=3;
571                    end
572                    ctable(ivec)=corrmax/sum_square;% correlation value
573                catch ME
574                    F(ivec)=3;
575                end
576            else
577                F(ivec)=3;
578            end
579        end
580    end
581   
582    %Create the vector matrix x, y, u, v
583end
584result_conv=result_conv*corrmax/(255*sum_square);% keep the last correlation matrix for output
585
586%------------------------------------------------------------------------
587% --- Find the maximum of the correlation function after interpolation
588function [vector,F] = SUBPIXGAUSS (result_conv,x,y)
589%------------------------------------------------------------------------
590vector=[0 0]; %default
591F=0;
592[npy,npx]=size(result_conv);
593
594% if (x <= (size(result_conv,1)-1)) && (y <= (size(result_conv,1)-1)) && (x >= 1) && (y >= 1)
595    %the following 8 lines are copyright (c) 1998, Uri Shavit, Roi Gurka, Alex Liberzon, Technion ï¿œ Israel Institute of Technology
596    %http://urapiv.wordpress.com
597    peaky = y;
598    if y <= npy-1 && y >= 1
599        f0 = log(result_conv(y,x));
600        f1 = real(log(result_conv(y-1,x)));
601        f2 = real(log(result_conv(y+1,x)));
602        peaky = peaky+ (f1-f2)/(2*f1-4*f0+2*f2);
603    else
604        F=-2; % warning flag for vector truncated by the limited search box
605    end
606    peakx=x;
607    if x <= npx-1 && x >= 1
608        f0 = log(result_conv(y,x));
609        f1 = real(log(result_conv(y,x-1)));
610        f2 = real(log(result_conv(y,x+1)));
611        peakx = peakx+ (f1-f2)/(2*f1-4*f0+2*f2);
612    else
613        F=-2; % warning flag for vector truncated by the limited search box
614    end
615    vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
616
617%------------------------------------------------------------------------
618% --- Find the maximum of the correlation function after interpolation
619function [vector,F] = SUBPIX2DGAUSS (result_conv,x,y)
620%------------------------------------------------------------------------
621vector=[0 0]; %default
622F=-2;
623peaky=y;
624peakx=x;
625[npy,npx]=size(result_conv);
626if (x <= npx-1) && (y <= npy-1) && (x >= 1) && (y >= 1)
627    F=0;
628    for i=-1:1
629        for j=-1:1
630            %following 15 lines based on
631            %H. Nobach ï¿œ M. Honkanen (2005)
632            %Two-dimensional Gaussian regression for sub-pixel displacement
633            %estimation in particle image velocimetry or particle position
634            %estimation in particle tracking velocimetry
635            %Experiments in Fluids (2005) 38: 511ï¿œ515
636            c10(j+2,i+2)=i*log(result_conv(y+j, x+i));
637            c01(j+2,i+2)=j*log(result_conv(y+j, x+i));
638            c11(j+2,i+2)=i*j*log(result_conv(y+j, x+i));
639            c20(j+2,i+2)=(3*i^2-2)*log(result_conv(y+j, x+i));
640            c02(j+2,i+2)=(3*j^2-2)*log(result_conv(y+j, x+i));
641        end
642    end
643    c10=(1/6)*sum(sum(c10));
644    c01=(1/6)*sum(sum(c01));
645    c11=(1/4)*sum(sum(c11));
646    c20=(1/6)*sum(sum(c20));
647    c02=(1/6)*sum(sum(c02));
648    deltax=(c11*c01-2*c10*c02)/(4*c20*c02-c11^2);
649    deltay=(c11*c10-2*c01*c20)/(4*c20*c02-c11^2);
650    if abs(deltax)<1
651        peakx=x+deltax;
652    end
653    if abs(deltay)<1
654        peaky=y+deltay;
655    end
656end
657vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
658
659%'RUN_FIX': function for fixing velocity fields:
660%-----------------------------------------------
661% RUN_FIX(filename,field,flagindex,thresh_vecC,thresh_vel,iter,flag_mask,maskname,fileref,fieldref)
662%
663%filename: name of the netcdf file (used as input and output)
664%field: structure specifying the names of the fields to fix (depending on civ1 or civ2)
665    %.vel_type='civ1' or 'civ2';
666    %.nb=name of the dimension common to the field to fix ('nb_vectors' for civ1);
667    %.fixflag=name of fix flag variable ('vec_FixFlag' for civ1)
668%flagindex: flag specifying which values of vec_f are removed:
669        % if flagindex(1)=1: vec_f=-2 vectors are removed
670        % if flagindex(2)=1: vec_f=3 vectors are removed
671        % if flagindex(3)=1: vec_f=2 vectors are removed (if iter=1) or vec_f=4 vectors are removed (if iter=2)
672%iter=1 for civ1 fields and iter=2 for civ2 fields
673%thresh_vecC: threshold in the image correlation vec_C
674%flag_mask: =1 mask used to remove vectors (0 else)
675%maskname: name of the mask image file for fix
676%thresh_vel: threshold on velocity, or on the difference with the reference file fileref if exists
677%inf_sup=1: remove values smaller than threshold thresh_vel, =2, larger than threshold
678%fileref: .nc file name for a reference velocity (='': refrence 0 used)
679%fieldref: 'civ1','filter1'...feld used in fileref
680
681function FF=fix(Param,F,C,U,V,X,Y)
682FF=zeros(size(F));%default
683
684%criterium on warn flags
685FlagName={'CheckFmin2','CheckF2','CheckF3','CheckF4'};
686FlagVal=[-2 2 3 4];
687for iflag=1:numel(FlagName)
688    if isfield(Param,FlagName{iflag}) && Param.(FlagName{iflag})
689        FF=(FF==1| F==FlagVal(iflag));
690    end
691end
692%criterium on correlation values
693if isfield (Param,'MinCorr')
694    FF=FF==1 | C<Param.MinCorr;
695end
696if (isfield(Param,'MinVel')&&~isempty(Param.MinVel))||(isfield (Param,'MaxVel')&&~isempty(Param.MaxVel))
697    Umod= U.*U+V.*V;
698    if isfield (Param,'MinVel')&&~isempty(Param.MinVel)
699        FF=FF==1 | Umod<(Param.MinVel*Param.MinVel);
700    end
701    if isfield (Param,'MaxVel')&&~isempty(Param.MaxVel)
702        FF=FF==1 | Umod>(Param.MaxVel*Param.MaxVel);
703    end
704end
705
706
707
708
709
710
711
Note: See TracBrowser for help on using the repository browser.