source: trunk/src/civ_matlab.m @ 456

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

bugs repaired

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