source: trunk/src/civ_matlab.m @ 690

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

various bugs corrected. Introduction of campaign lists in Open menu

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