source: trunk/src/civ_matlab.m @ 441

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

rationalisation of civ parameter names
civ results now put in a directory outside the image directory.
*xml file management needs to be updated

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