source: trunk/src/civ_matlab.m @ 386

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

several bugs corrected
for dealing with color movies

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