source: trunk/src/civ_matlab.m @ 406

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

bugs corrected in civ_matlab and object projection

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