source: trunk/src/series/civ_series.m @ 606

Last change on this file since 606 was 606, checked in by sommeria, 11 years ago

bug on compilation solved (still to test with transform_field fct). faster browser inrtroduced. Various bug corrections

File size: 42.8 KB
Line 
1%'civ_series': PIV function activated by the general GUI series
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_series(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_series(Param,ncfile)
41errormsg='';
42path_series=fileparts(which('series'));
43addpath(fullfile(path_series,'series'))
44%% set the input elements needed on the GUI series when the action is selected in the menu ActionName or InputTable refreshed
45if isstruct(Param) && isequal(Param.Action.RUN,0)
46    Data=civ_input(Param);% introduce the civ parameters using the GUI civ_input
47    Data.Program=mfilename;%gives the name of the current function
48    Data.AllowInputSort='off';% allow alphabetic sorting of the list of input file SubDir (options 'off'/'on', 'off' by default)
49    Data.WholeIndexRange='off';% prescribes the file index ranges from min to max (options 'off'/'on', 'off' by default)
50    Data.NbSlice='off'; %nbre of slices ('off' by default)
51    Data.VelType='off';% menu for selecting the velocity type (options 'off'/'one'/'two',  'off' by default)
52    Data.FieldName='off';% menu for selecting the field (s) in the input file(options 'off'/'one'/'two', 'off' by default)
53    Data.FieldTransform = 'off';%can use a transform function
54    Data.ProjObject='off';%can use projection object(option 'off'/'on',
55    Data.Mask='off';%can use mask option   (option 'off'/'on', 'off' by default)
56    Data.OutputDirExt='.civ';%set the output dir extension
57    filecell=get_file_series(Param);%check existence of the first input file
58    if ~exist(filecell{1,1},'file')
59        msgbox_uvmat('WARNING','the first input file does not exist')
60    end
61    return
62end
63
64%% read input parameters from an xml file if input is a file name (batch mode)
65checkrun=1;
66if ischar(Param)
67    Param=xml2struct(Param);% read Param as input file (batch case)
68    checkrun=0;
69end
70
71%% input files and indexing
72NbField=1;
73if isfield(Param,'InputTable')
74    RootPath=Param.InputTable{1,1};
75    RootFile=Param.InputTable{1,3};
76    SubDir=Param.InputTable{1,2};
77    NomType=Param.InputTable{1,4};
78    FileExt=Param.InputTable{1,5};
79    PairCiv1=Param.ActionInput.PairIndices.ListPairCiv1;
80    PairCiv2='';
81    if isfield(Param.ActionInput.PairIndices,'ListPairCiv2')
82        PairCiv2=Param.ActionInput.PairIndices.ListPairCiv2;
83    end
84    MaxIndex=cell2mat(Param.IndexRange.MaxIndex);
85    MinIndex=cell2mat(Param.IndexRange.MinIndex);
86    [filecell,i_series,tild,j_series]=get_file_series(Param);
87    [i1_series_Civ1,i2_series_Civ1,j1_series_Civ1,j2_series_Civ1,check_bounds,NomTypeNc]=...
88        find_pair_indices(PairCiv1,i_series{1},j_series{1},MinIndex,MaxIndex);
89    if ~isempty(PairCiv2)
90        [i1_series_Civ2,i2_series_Civ2,j1_series_Civ2,j2_series_Civ2,check_bounds_Civ2]=...
91            find_pair_indices(PairCiv2,i_series{1},j_series{1},MinIndex,MaxIndex);
92        check_bounds=check_bounds | check_bounds_Civ2;
93    end
94    i1_series_Civ1=i1_series_Civ1(~check_bounds);
95    i2_series_Civ1=i2_series_Civ1(~check_bounds);
96    j1_series_Civ1=j1_series_Civ1(~check_bounds);
97    j2_series_Civ1=j2_series_Civ1(~check_bounds);
98    if ~isempty(j1_series_Civ1)
99        FrameIndex_A_Civ1=j1_series_Civ1;
100        FrameIndex_B_Civ1=j2_series_Civ1;
101    else
102        FrameIndex_A_Civ1=i1_series_Civ1;
103        FrameIndex_B_Civ1=i2_series_Civ1;
104    end
105    if ~isempty(PairCiv2)
106        i1_series_Civ2=i1_series_Civ2(~check_bounds);
107        i2_series_Civ2=i2_series_Civ2(~check_bounds);
108        j1_series_Civ2=j1_series_Civ2(~check_bounds);
109        j2_series_Civ2=j2_series_Civ2(~check_bounds);
110        if ~isempty(j1_series_Civ2)
111            FrameIndex_A_Civ2=j1_series_Civ2;
112            FrameIndex_B_Civ2=j2_series_Civ2;
113        else
114            FrameIndex_A_Civ2=i1_series_Civ2;
115            FrameIndex_B_Civ2=i2_series_Civ2;
116        end
117    end
118   
119    NbField=numel(i1_series_Civ1);
120    [FileType_A,FileInfo,MovieObject_A]=get_file_type(filecell{1,1});
121    FileType_B=FileType_A;
122    MovieObject_B=MovieObject_A;
123    if size(filecell,1)>=2 && ~strcmp(filecell{1,1},filecell{2,1})
124        [FileType_B,FileInfo,MovieObject_B]=get_file_type(filecell{2,1});
125    end
126end
127
128
129%% Output directory
130OutputDir=[Param.OutputSubDir Param.OutputDirExt];
131
132Data.ListGlobalAttribute={'Conventions','Program','CivStage'};
133Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
134Data.Program='civ_series';
135Data.CivStage=0;%default
136ListVarCiv1={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F'}; %variables to read
137ListVarFix1={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_F','Civ1_FF'};
138mask='';
139maskname='';%default
140check_civx=0;%default
141check_civ1=0;%default
142check_patch1=0;%default
143
144%% get timing from the ImaDoc file or input video
145[XmlData,NbSlice_calib,time,errormsg]=read_multimadoc(RootPath,SubDir,RootFile,FileExt,i1_series,i2_series,j1_series,j2_series);
146%TODO: get time_A and time_B
147% case of movies TODO TODO TODO
148if isempty(time) && (strcmp(FileType,'video') || strcmp(FileType,'mmreader'))
149    set(handles.ListPairMode,'Value',1);
150    dt=1/get(MovieObject,'FrameRate');%time interval between successive frames
151    if strcmp(NomTypeIma,'*')
152        set(handles.ListPairMode,'String',{'series(Di)'})
153        MaxIndex_i=get(MovieObject,'NumberOfFrames');
154        time=(dt*(0:MaxIndex_i-1))';%list of image times
155    else
156        set(handles.ListPairMode,'String',[{'series(Dj)'};{'series(Di)'}])
157        MaxIndex_i=max(i1_series(i1_series>0));
158        MaxIndex_j=get(MovieObject,'NumberOfFrames');
159        time=ones(MaxIndex_i,1)*(dt*(0:MaxIndex_j-1));%list of image times
160        enable_j(handles,'on')
161    end
162    TimeUnit='s';
163%%%%% MAIN LOOP %%%%%%
164
165MovieObject_A=[];
166for ifield=1:NbField
167   
168    %% Civ1
169    if isfield (Param.ActionInput,'Civ1')
170        par_civ1=Param.ActionInput.Civ1;
171        if isfield(par_civ1,'reverse_pair')% A REVOIR
172            if par_civ1.reverse_pair
173                if ischar(par_civ1.ImageB)
174                    temp=par_civ1.ImageA;
175                    par_civ1.ImageA=imread(par_civ1.ImageB);
176                end
177                if ischar(temp)
178                    par_civ1.ImageB=imread(temp);
179                end
180            end
181        else
182            ImageName_A=fullfile_uvmat(RootPath,SubDir,RootFile,FileExt,NomType,i1_series_Civ1(ifield),[],j1_series_Civ1(ifield));
183            [par_civ1.ImageA,MovieObject_A] = read_image(ImageName_A,FileType_A,MovieObject_A,FrameIndex_A_Civ1(ifield));
184            ImageName_B=fullfile_uvmat(RootPath,SubDir,RootFile,FileExt,NomType,i2_series_Civ1(ifield),[],j2_series_Civ1(ifield));
185            [par_civ1.ImageB,MovieObject_B] = read_image(ImageName_B,FileType_B,MovieObject_B,FrameIndex_B_Civ1(ifield));
186        end
187        ncfile=fullfile_uvmat(RootPath,OutputDir,RootFile,'.nc',NomTypeNc,i1_series_Civ1(ifield),i2_series_Civ1(ifield),...
188            j1_series_Civ1(ifield),j2_series_Civ1(ifield));
189        par_civ1.ImageWidth=FileInfo.Width;
190        par_civ1.ImageHeight=FileInfo.Height;
191        list_param=(fieldnames(Param.ActionInput.Civ1))';
192        Civ1_param=regexprep(list_param,'^.+','Civ1_$0');% insert 'Civ1_' before  each string in list_param
193        Civ1_param=[{'Civ1_ImageA','Civ1_ImageB','Civ1_Time','Civ1_Dt'} Civ1_param]; %insert the names of the two input images
194        %indicate the values of all the global attributes in the output data
195        Data.Civ1_ImageA=ImageName_A;
196        Data.Civ1_ImageB=ImageName_B;
197        Data.Civ1_Time=((time(i2_civ1(ifile)+1,j2_civ1(j)+1)+time(i1_civ1(ifile)+1,j1_civ1(j)+1))/2);
198        Data.Civ1_Dt=(time(i2_civ1(ifile)+1,j2_civ1(j)+1)-time(i1_civ1(ifile)+1,j1_civ1(j)+1));
199        for ilist=1:length(list_param)
200            Data.(Civ1_param{4+ilist})=Param.ActionInput.Civ1.(list_param{ilist});
201        end
202        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ1_param];
203        Data.CivStage=1;
204       
205        % set the list of variables
206        Data.ListVarName={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_F','Civ1_C'};%  cell array containing the names of the fields to record
207        Data.VarDimName={'nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1'};
208        Data.VarAttribute{1}.Role='coord_x';
209        Data.VarAttribute{2}.Role='coord_y';
210        Data.VarAttribute{3}.Role='vector_x';
211        Data.VarAttribute{4}.Role='vector_y';
212        Data.VarAttribute{5}.Role='warnflag';
213       
214        if strcmp(Param.ActionInput.ListCompareMode, 'PIV volume')
215            Data.ListVarName=[Data.ListVarName 'Civ1_Z'];
216            Data.Civ1_X=[];Data.Civ1_Y=[];Data.Civ1_Z=[];
217            Data.Civ1_U=[];Data.Civ1_V=[];Data.Civ1_C=[];Data.Civ1_F=[];
218            for ivol=1:NbSlice
219                % caluclate velocity data (y and v in indices, reverse to y component)
220                [xtable ytable utable vtable ctable F result_conv errormsg] = civ (par_civ1);
221                if ~isempty(errormsg)
222                    return
223                end
224                Data.Civ1_X=[Data.Civ1_X reshape(xtable,[],1)];
225                Data.Civ1_Y=[Data.Civ1_Y reshape(Param.Civ1.ImageHeight-ytable+1,[],1)];
226                Data.Civ1_Z=[Data.Civ1_Z ivol*ones(numel(xtable),1)];% z=image index in image coordinates
227                Data.Civ1_U=[Data.Civ1_U reshape(utable,[],1)];
228                Data.Civ1_V=[Data.Civ1_V reshape(-vtable,[],1)];
229                Data.Civ1_C=[Data.Civ1_C reshape(ctable,[],1)];
230                Data.Civ1_F=[Data.Civ1_C reshape(F,[],1)];
231            end
232        else %usual PIV
233            % caluclate velocity data (y and v in indices, reverse to y component)
234            [xtable ytable utable vtable ctable F result_conv errormsg] = civ (par_civ1);
235            if ~isempty(errormsg)
236                return
237            end
238            Data.Civ1_X=reshape(xtable,[],1);
239            Data.Civ1_Y=reshape(par_civ1.ImageHeight-ytable+1,[],1);
240            Data.Civ1_U=reshape(utable,[],1);
241            Data.Civ1_V=reshape(-vtable,[],1);
242            Data.Civ1_C=reshape(ctable,[],1);
243            Data.Civ1_F=reshape(F,[],1);
244        end
245    else
246        if exist('ncfile','var')
247            CivFile=ncfile;
248        elseif isfield(Param.Patch1,'CivFile')
249            CivFile=Param.Patch1.CivFile;
250        end
251        Data=nc2struct(CivFile,'ListGlobalAttribute','absolut_time_T0'); %look for the constant 'absolut_time_T0' to detect old civx data format
252        if isfield(Data,'Txt')
253            errormsg=Data.Txt;
254            return
255        end
256        if ~isempty(Data.absolut_time_T0')%read civx file
257            check_civx=1;% test for old civx data format
258            [Data,vardetect,ichoice]=nc2struct(CivFile);%read the variables in the netcdf file
259        else
260            Data=nc2struct(CivFile);%read civ1 and fix1 data in the existing netcdf file
261        end
262    end
263   
264    %% Fix1
265    if isfield (Param.ActionInput,'Fix1')
266        ListFixParam=fieldnames(Param.ActionInput.Fix1);
267        for ilist=1:length(ListFixParam)
268            ParamName=ListFixParam{ilist};
269            ListName=['Fix1_' ParamName];
270            eval(['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];'])
271            eval(['Data.' ListName '=Param.ActionInput.Fix1.' ParamName ';'])
272        end
273        if check_civx
274            if ~isfield(Data,'fix')
275                Data.ListGlobalAttribute=[Data.ListGlobalAttribute 'fix'];
276                Data.fix=1;
277                Data.ListVarName=[Data.ListVarName {'vec_FixFlag'}];
278                Data.VarDimName=[Data.VarDimName {'nb_vectors'}];
279            end
280            Data.vec_FixFlag=fix(Param.ActionInput.Fix1,Data.vec_F,Data.vec_C,Data.vec_U,Data.vec_V,Data.vec_X,Data.vec_Y);
281        else
282            Data.ListVarName=[Data.ListVarName {'Civ1_FF'}];
283            Data.VarDimName=[Data.VarDimName {'nb_vec_1'}];
284            nbvar=length(Data.ListVarName);
285            Data.VarAttribute{nbvar}.Role='errorflag';
286            Data.Civ1_FF=fix(Param.ActionInput.Fix1,Data.Civ1_F,Data.Civ1_C,Data.Civ1_U,Data.Civ1_V);
287            Data.CivStage=2;
288        end
289    end
290    %% Patch1
291    if isfield (Param.ActionInput,'Patch1')
292        if check_civx
293            errormsg='Civ Matlab input needed for patch';
294            return
295        end
296       
297        Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Patch1_Rho','Patch1_Threshold','Patch1_SubDomain'}];
298        Data.Patch1_FieldSmooth=Param.ActionInput.Patch1.FieldSmooth;
299        Data.Patch1_MaxDiff=Param.ActionInput.Patch1.MaxDiff;
300        Data.Patch1_SubDomainSize=Param.ActionInput.Patch1.SubDomainSize;
301        nbvar=length(Data.ListVarName);
302        Data.ListVarName=[Data.ListVarName {'Civ1_U_smooth','Civ1_V_smooth','Civ1_SubRange','Civ1_NbCentres','Civ1_Coord_tps','Civ1_U_tps','Civ1_V_tps'}];
303        Data.VarDimName=[Data.VarDimName {'nb_vec_1','nb_vec_1',{'nb_coord','nb_bounds','nb_subdomain_1'},'nb_subdomain_1',...
304            {'nb_tps_1','nb_coord','nb_subdomain_1'},{'nb_tps_1','nb_subdomain_1'},{'nb_tps_1','nb_subdomain_1'}}];
305        Data.VarAttribute{nbvar+1}.Role='vector_x';
306        Data.VarAttribute{nbvar+2}.Role='vector_y';
307        Data.VarAttribute{nbvar+5}.Role='coord_tps';
308        Data.VarAttribute{nbvar+6}.Role='vector_x';
309        Data.VarAttribute{nbvar+7}.Role='vector_y';
310        Data.Civ1_U_smooth=zeros(size(Data.Civ1_X));
311        Data.Civ1_V_smooth=zeros(size(Data.Civ1_X));
312        if isfield(Data,'Civ1_FF')
313            ind_good=find(Data.Civ1_FF==0);
314        else
315            ind_good=1:numel(Data.Civ1_X);
316        end
317        [Data.Civ1_SubRange,Data.Civ1_NbCentres,Data.Civ1_Coord_tps,Data.Civ1_U_tps,Data.Civ1_V_tps,tild,Ures, Vres,tild,FFres]=...
318            filter_tps([Data.Civ1_X(ind_good) Data.Civ1_Y(ind_good)],Data.Civ1_U(ind_good),Data.Civ1_V(ind_good),[],Data.Patch1_SubDomainSize,Data.Patch1_FieldSmooth,Data.Patch1_MaxDiff);
319        Data.Civ1_U_smooth(ind_good)=Ures;
320        Data.Civ1_V_smooth(ind_good)=Vres;
321        Data.Civ1_FF(ind_good)=FFres;
322        Data.CivStage=3;
323    end
324   
325    %% Civ2
326    if isfield (Param.ActionInput,'Civ2')
327        par_civ2=Param.ActionInput.Civ2;
328        par_civ2.ImageA=[];
329        par_civ2.ImageB=[];
330        %         if ~isfield(Param.Civ1,'ImageA')
331        ImageName_A_Civ2=fullfile_uvmat(RootPath,SubDir,RootFile,FileExt,NomType,i1_series_Civ2(ifield),[],j1_series_Civ2(ifield));
332
333        if strcmp(ImageName_A_Civ2,ImageName_A) && isequal(FrameIndex_A_Civ1(ifield),FrameIndex_A_Civ2)
334            par_civ2.ImageA=par_civ1.ImageA;
335        else
336            [par_civ2.ImageA,MovieObject_A] = read_image(ImageName_A,FileType_A,MovieObject_A,FrameIndex_A_Civ2);
337        end
338        ImageName_B_Civ2=fullfile_uvmat(RootPath,SubDir,RootFile,FileExt,NomType,i2_series_Civ2(ifield),[],j2_series_Civ2(ifield));
339        if strcmp(ImageName_B_Civ2,ImageName_B) && isequal(FrameIndex_B_Civ1(ifield),FrameIndex_B_Civ2)
340            par_civ2.ImageB=par_civ1.ImageB;
341        else
342            [par_civ2.ImageB,MovieObject_B] = read_image(ImageName_B,FileType_B,MovieObject_B,FrameIndex_B_Civ2);
343        end     
344       
345        ncfile=fullfile_uvmat(RootPath,OutputDir,RootFile,'.nc',NomTypeNc,i1_series_Civ2(ifield),i2_series_Civ2(ifield),...
346            j1_series_Civ2(ifield),j2_series_Civ2(ifield));
347        par_civ2.ImageWidth=FileInfo.Width;
348        par_civ2.ImageHeight=FileInfo.Height;
349       
350        if isfield(par_civ2,'Grid')% grid points set as input file
351            if ischar(par_civ2.Grid)%read the grid file if the input is a file name
352                par_civ2.Grid=dlmread(par_civ2.Grid);
353                par_civ2.Grid(1,:)=[];%the first line must be removed (heading in the grid file)
354            end
355        else% automatic grid
356            minix=floor(par_civ2.Dx/2)-0.5;
357            maxix=minix+par_civ2.Dx*floor((par_civ2.ImageWidth-1)/par_civ2.Dx);
358            miniy=floor(par_civ2.Dy/2)-0.5;
359            maxiy=minix+par_civ2.Dy*floor((par_civ2.ImageHeight-1)/par_civ2.Dy);
360            [GridX,GridY]=meshgrid(minix:par_civ2.Dx:maxix,miniy:par_civ2.Dy:maxiy);
361            par_civ2.Grid(:,1)=reshape(GridX,[],1);
362            par_civ2.Grid(:,2)=reshape(GridY,[],1);
363        end
364        Shiftx=zeros(size(par_civ2.Grid,1),1);% shift expected from civ1 data
365        Shifty=zeros(size(par_civ2.Grid,1),1);
366        nbval=zeros(size(par_civ2.Grid,1),1);
367        if par_civ2.CheckDeformation
368            DUDX=zeros(size(par_civ2.Grid,1),1);
369            DUDY=zeros(size(par_civ2.Grid,1),1);
370            DVDX=zeros(size(par_civ2.Grid,1),1);
371            DVDY=zeros(size(par_civ2.Grid,1),1);
372        end
373        NbSubDomain=size(Data.Civ1_SubRange,3);
374        % get the guess from patch1
375        for isub=1:NbSubDomain
376            nbvec_sub=Data.Civ1_NbCentres(isub);
377            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));
378            epoints = [GridX(ind_sel) GridY(ind_sel)];% coordinates of interpolation sites
379            ctrs=Data.Civ1_Coord_tps(1:nbvec_sub,:,isub) ;%(=initial points) ctrs
380            nbval(ind_sel)=nbval(ind_sel)+1;% records the number of values for eacn interpolation point (in case of subdomain overlap)
381            EM = tps_eval(epoints,ctrs);
382            Shiftx(ind_sel)=Shiftx(ind_sel)+EM*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
383            Shifty(ind_sel)=Shifty(ind_sel)+EM*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
384            if par_civ2.CheckDeformation
385                [EMDX,EMDY] = tps_eval_dxy(epoints,ctrs);%2D matrix of distances between extrapolation points epoints and spline centres (=site points) ctrs
386                DUDX(ind_sel)=DUDX(ind_sel)+EMDX*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
387                DUDY(ind_sel)=DUDY(ind_sel)+EMDY*Data.Civ1_U_tps(1:nbvec_sub+3,isub);
388                DVDX(ind_sel)=DVDX(ind_sel)+EMDX*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
389                DVDY(ind_sel)=DVDY(ind_sel)+EMDY*Data.Civ1_V_tps(1:nbvec_sub+3,isub);
390            end
391        end
392        mask='';
393        if par_civ2.CheckMask&&~isempty(par_civ2.Mask)&& ~strcmp(maskname,par_civ2.Mask)% mask exist, not already read in civ1
394            mask=imread(par_civ2.Mask);
395        end
396        ibx2=ceil(par_civ2.CorrBoxSize(1)/2);
397        iby2=ceil(par_civ2.CorrBoxSize(2)/2);
398        par_civ2.SearchBoxSize(1)=2*ibx2+9;% search ara +-4 pixels around the guess
399        par_civ2.SearchBoxSize(2)=2*iby2+9;
400        par_civ2.SearchBoxShift=[Shiftx(nbval>=1)./nbval(nbval>=1) Shifty(nbval>=1)./nbval(nbval>=1)];
401        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
402        if par_civ2.CheckDeformation
403            par_civ2.DUDX=DUDX./nbval;
404            par_civ2.DUDY=DUDY./nbval;
405            par_civ2.DVDX=DVDX./nbval;
406            par_civ2.DVDY=DVDY./nbval;
407        end
408        % caluclate velocity data (y and v in indices, reverse to y component)
409        [xtable ytable utable vtable ctable F] = civ (par_civ2);
410
411        list_param=(fieldnames(Param.ActionInput.Civ2))';
412        Civ2_param=regexprep(list_param,'^.+','Civ2_$0');% insert 'Civ2_' before  each string in list_param
413        Civ2_param=[{'Civ2_ImageA','Civ2_ImageB','Civ2_Time','Civ2_Dt'} Civ2_param]; %insert the names of the two input images
414        %indicate the values of all the global attributes in the output data
415        Data.Civ2_ImageA=ImageName_A;
416        Data.Civ2_ImageB=ImageName_B;
417        Data.Civ2_Time=1;
418        Data.Civ2_Dt=1;
419        for ilist=1:length(list_param)
420            Data.(Civ2_param{4+ilist})=Param.ActionInput.Civ2.(list_param{ilist});
421        end
422        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ2_param];
423       
424        nbvar=numel(Data.ListVarName);
425        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
426        Data.VarDimName=[Data.VarDimName {'nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2'}];
427        Data.VarAttribute{nbvar+1}.Role='coord_x';
428        Data.VarAttribute{nbvar+2}.Role='coord_y';
429        Data.VarAttribute{nbvar+3}.Role='vector_x';
430        Data.VarAttribute{nbvar+4}.Role='vector_y';
431        Data.VarAttribute{nbvar+5}.Role='warnflag';
432        Data.Civ2_X=reshape(xtable,[],1);
433        Data.Civ2_Y=reshape(size(par_civ2.ImageA,1)-ytable+1,[],1);
434        Data.Civ2_U=reshape(utable,[],1);
435        Data.Civ2_V=reshape(-vtable,[],1);
436        Data.Civ2_C=reshape(ctable,[],1);
437        Data.Civ2_F=reshape(F,[],1);
438        Data.CivStage=Data.CivStage+1;
439    end
440   
441    %% Fix2
442    if isfield (Param.ActionInput,'Fix2')
443        ListFixParam=fieldnames(Param.ActionInput.Fix2);
444        for ilist=1:length(ListFixParam)
445            ParamName=ListFixParam{ilist};
446            ListName=['Fix2_' ParamName];
447            eval(['Data.ListGlobalAttribute=[Data.ListGlobalAttribute ''' ParamName '''];'])
448            eval(['Data.' ListName '=Param.ActionInput.Fix2.' ParamName ';'])
449        end
450        if check_civx
451            if ~isfield(Data,'fix2')
452                Data.ListGlobalAttribute=[Data.ListGlobalAttribute 'fix2'];
453                Data.fix2=1;
454                Data.ListVarName=[Data.ListVarName {'vec2_FixFlag'}];
455                Data.VarDimName=[Data.VarDimName {'nb_vectors2'}];
456            end
457            Data.vec_FixFlag=fix(Param.Fix2,Data.vec2_F,Data.vec2_C,Data.vec2_U,Data.vec2_V,Data.vec2_X,Data.vec2_Y);
458        else
459            Data.ListVarName=[Data.ListVarName {'Civ2_FF'}];
460            Data.VarDimName=[Data.VarDimName {'nb_vec_2'}];
461            nbvar=length(Data.ListVarName);
462            Data.VarAttribute{nbvar}.Role='errorflag';
463            Data.Civ2_FF=fix(Param.ActionInput.Fix2,Data.Civ2_F,Data.Civ2_C,Data.Civ2_U,Data.Civ2_V);
464            Data.CivStage=Data.CivStage+1;
465        end
466       
467    end
468   
469    %% Patch2
470    if isfield (Param.ActionInput,'Patch2')
471        Data.ListGlobalAttribute=[Data.ListGlobalAttribute {'Patch2_Rho','Patch2_Threshold','Patch2_SubDomain'}];
472        Data.Patch2_FieldSmooth=Param.ActionInput.Patch2.FieldSmooth;
473        Data.Patch2_MaxDiff=Param.ActionInput.Patch2.MaxDiff;
474        Data.Patch2_SubDomainSize=Param.ActionInput.Patch2.SubDomainSize;
475        nbvar=length(Data.ListVarName);
476        Data.ListVarName=[Data.ListVarName {'Civ2_U_smooth','Civ2_V_smooth','Civ2_SubRange','Civ2_NbCentres','Civ2_Coord_tps','Civ2_U_tps','Civ2_V_tps'}];
477        Data.VarDimName=[Data.VarDimName {'nb_vec_2','nb_vec_2',{'nb_coord','nb_bounds','nb_subdomain_2'},{'nb_subdomain_2'},...
478            {'nb_tps_2','nb_coord','nb_subdomain_2'},{'nb_tps_2','nb_subdomain_2'},{'nb_tps_2','nb_subdomain_2'}}];
479       
480        Data.VarAttribute{nbvar+1}.Role='vector_x';
481        Data.VarAttribute{nbvar+2}.Role='vector_y';
482        Data.VarAttribute{nbvar+5}.Role='coord_tps';
483        Data.VarAttribute{nbvar+6}.Role='vector_x';
484        Data.VarAttribute{nbvar+7}.Role='vector_y';
485        Data.Civ2_U_smooth=zeros(size(Data.Civ2_X));
486        Data.Civ2_V_smooth=zeros(size(Data.Civ2_X));
487        if isfield(Data,'Civ2_FF')
488            ind_good=find(Data.Civ2_FF==0);
489        else
490            ind_good=1:numel(Data.Civ2_X);
491        end
492        [Data.Civ2_SubRange,Data.Civ2_NbCentres,Data.Civ2_Coord_tps,Data.Civ2_U_tps,Data.Civ2_V_tps,tild,Ures, Vres,tild,FFres]=...
493            filter_tps([Data.Civ2_X(ind_good) Data.Civ2_Y(ind_good)],Data.Civ2_U(ind_good),Data.Civ2_V(ind_good),[],Data.Patch2_SubDomainSize,Data.Patch2_FieldSmooth,Data.Patch2_MaxDiff);
494        Data.Civ2_U_smooth(ind_good)=Ures;
495        Data.Civ2_V_smooth(ind_good)=Vres;
496        Data.Civ2_FF(ind_good)=FFres;
497        Data.CivStage=Data.CivStage+1;
498    end
499   
500    %% write result in a netcdf file if requested
501    if exist('ncfile','var')
502        errormsg=struct2nc(ncfile,Data);
503        if isempty(errormsg)
504            disp([ncfile ' written'])
505        else
506            disp(errormsg)
507        end
508    end
509end
510
511% 'civ': function piv.m adapted from PIVlab http://pivlab.blogspot.com/
512%--------------------------------------------------------------------------
513% function [xtable ytable utable vtable typevector] = civ (image1,image2,ibx,iby step, subpixfinder, mask, roi)
514%
515% OUTPUT:
516% xtable: set of x coordinates
517% ytable: set of y coordiantes
518% utable: set of u displacements (along x)
519% vtable: set of v displacements (along y)
520% ctable: max image correlation for each vector
521% typevector: set of flags, =1 for good, =0 for NaN vectors
522%
523%INPUT:
524% par_civ: structure of input parameters, with fields:
525%  .CorrBoxSize
526%  .SearchBoxSize
527%  .SearchBoxShift
528%  .ImageHeight
529%  .ImageWidth
530%  .Dx, Dy
531%  .Grid
532%  .Mask
533%  .MinIma
534%  .MaxIma
535%  .image1:first image (matrix)
536% image2: second image (matrix)
537% ibx2,iby2: half size of the correlation box along x and y, in px (size=(2*iby2+1,2*ibx2+1)
538% isx2,isy2: half size of the search box along x and y, in px (size=(2*isy2+1,2*isx2+1)
539% shiftx, shifty: shift of the search box (in pixel index, yshift reversed)
540% step: mesh of the measurement points (in px)
541% subpixfinder=1 or 2 controls the curve fitting of the image correlation
542% mask: =[] for no mask
543% roi: 4 element vector defining a region of interest: x position, y position, width, height, (in image indices), for the whole image, roi=[];
544function [xtable ytable utable vtable ctable F result_conv errormsg] = civ (par_civ)
545%this funtion performs the DCC PIV analysis. Recent window-deformation
546%methods perform better and will maybe be implemented in the future.
547
548%% prepare measurement grid
549if isfield(par_civ,'Grid')% grid points set as input
550    if ischar(par_civ.Grid)%read the drid file if the input is a file name
551        par_civ.Grid=dlmread(par_civ.Grid);
552        par_civ.Grid(1,:)=[];%the first line must be removed (heading in the grid file)
553    end
554else% automatic grid
555    minix=floor(par_civ.Dx/2)-0.5;
556    maxix=minix+par_civ.Dx*floor((par_civ.ImageWidth-1)/par_civ.Dx);
557    miniy=floor(par_civ.Dy/2)-0.5;
558    maxiy=minix+par_civ.Dy*floor((par_civ.ImageHeight-1)/par_civ.Dy);
559    [GridX,GridY]=meshgrid(minix:par_civ.Dx:maxix,miniy:par_civ.Dy:maxiy);
560    par_civ.Grid(:,1)=reshape(GridX,[],1);
561    par_civ.Grid(:,2)=reshape(GridY,[],1);
562end
563nbvec=size(par_civ.Grid,1);
564
565%% prepare correlation and search boxes
566ibx2=ceil(par_civ.CorrBoxSize(1)/2);
567iby2=ceil(par_civ.CorrBoxSize(2)/2);
568isx2=ceil(par_civ.SearchBoxSize(1)/2);
569isy2=ceil(par_civ.SearchBoxSize(2)/2);
570shiftx=round(par_civ.SearchBoxShift(:,1));
571shifty=-round(par_civ.SearchBoxShift(:,2));% sign minus because image j index increases when y decreases
572if numel(shiftx)==1% case of a unique shift for the whole field( civ1)
573    shiftx=shiftx*ones(nbvec,1);
574    shifty=shifty*ones(nbvec,1);
575end
576
577%% Default output
578xtable=par_civ.Grid(:,1);
579ytable=par_civ.Grid(:,2);
580utable=zeros(nbvec,1);
581vtable=zeros(nbvec,1);
582ctable=zeros(nbvec,1);
583F=zeros(nbvec,1);
584result_conv=[];
585errormsg='';
586
587%% prepare mask
588if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
589    if strcmp(par_civ.Mask,'all')
590        return    % get the grid only, no civ calculation
591    elseif ischar(par_civ.Mask)
592        par_civ.Mask=imread(par_civ.Mask);
593    end
594end
595check_MinIma=isfield(par_civ,'MinIma');% test for image luminosity threshold
596check_MaxIma=isfield(par_civ,'MaxIma') && ~isempty(par_civ.MaxIma);
597
598% %% prepare images
599% if isfield(par_civ,'reverse_pair')
600%     if par_civ.reverse_pair
601%         if ischar(par_civ.ImageB)
602%             temp=par_civ.ImageA;
603%             par_civ.ImageA=imread(par_civ.ImageB);
604%         end
605%         if ischar(temp)
606%             par_civ.ImageB=imread(temp);
607%         end
608%     end
609% else
610%     if ischar(par_civ.ImageA)
611%         par_civ.ImageA=imread(par_civ.ImageA);
612%     end
613%     if ischar(par_civ.ImageB)
614%         par_civ.ImageB=imread(par_civ.ImageB);
615%     end
616% end
617par_civ.ImageA=sum(double(par_civ.ImageA),3);%sum over rgb component for color images
618par_civ.ImageB=sum(double(par_civ.ImageB),3);
619[npy_ima npx_ima]=size(par_civ.ImageA);
620if ~isequal(size(par_civ.ImageB),[npy_ima npx_ima])
621    errormsg='image pair with unequal size';
622    return
623end
624
625%% Apply mask
626    % Convention for mask IDEAS TO IMPLEMENT ?
627    % mask >200 : velocity calculated
628    %  200 >=mask>150;velocity not calculated, interpolation allowed (bad spots)
629    % 150>=mask >100: velocity not calculated, nor interpolated
630    %  100>=mask> 20: velocity not calculated, impermeable (no flux through mask boundaries)
631    %  20>=mask: velocity=0
632checkmask=0;
633MinA=min(min(par_civ.ImageA));
634MinB=min(min(par_civ.ImageB));
635if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
636   checkmask=1;
637   if ~isequal(size(par_civ.Mask),[npy_ima npx_ima])
638        errormsg='mask must be an image with the same size as the images';
639        return
640   end
641  %  check_noflux=(par_civ.Mask<100) ;%TODO: to implement
642    check_undefined=(par_civ.Mask<200 & par_civ.Mask>=20 );
643    par_civ.ImageA(check_undefined)=MinA;% put image A to zero (i.e. the min image value) in the undefined  area
644    par_civ.ImageB(check_undefined)=MinB;% put image B to zero (i.e. the min image value) in the undefined  area
645end
646
647%% compute image correlations: MAINLOOP on velocity vectors
648corrmax=0;
649sum_square=1;% default
650mesh=1;% default
651CheckDecimal=isfield(par_civ,'CheckDecimal')&& par_civ.CheckDecimal==1;
652if CheckDecimal
653    mesh=0.2;%mesh in pixels for subpixel image interpolation
654    CheckDeformation=isfield(par_civ,'CheckDeformation')&& par_civ.CheckDeformation==1;
655end
656% vector=[0 0];%default
657for ivec=1:nbvec
658    iref=round(par_civ.Grid(ivec,1)+0.5);% xindex on the image A for the middle of the correlation box
659    jref=round(par_civ.ImageHeight-par_civ.Grid(ivec,2)+0.5);% yindex on the image B for the middle of the correlation box
660    %if ~(checkmask && par_civ.Mask(jref,iref)<=20) %velocity not set to zero by the black mask
661    %         if jref-iby2<1 || jref+iby2>par_civ.ImageHeight|| iref-ibx2<1 || iref+ibx2>par_civ.ImageWidth||...
662    %               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
663    %             F(ivec)=3;
664    %         else
665    F(ivec)=0;
666    subrange1_x=iref-ibx2:iref+ibx2;% x indices defining the first subimage
667    subrange1_y=jref-iby2:jref+iby2;% y indices defining the first subimage
668    subrange2_x=iref+shiftx(ivec)-isx2:iref+shiftx(ivec)+isx2;%x indices defining the second subimage
669    subrange2_y=jref+shifty(ivec)-isy2:jref+shifty(ivec)+isy2;%y indices defining the second subimage
670    image1_crop=MinA*ones(numel(subrange1_y),numel(subrange1_x));% default value=min of image A
671    image2_crop=MinA*ones(numel(subrange2_y),numel(subrange2_x));% default value=min of image A
672    check1_x=subrange1_x>=1 & subrange1_x<=par_civ.ImageWidth;% check which points in the subimage 1 are contained in the initial image 1
673    check1_y=subrange1_y>=1 & subrange1_y<=par_civ.ImageHeight;
674    check2_x=subrange2_x>=1 & subrange2_x<=par_civ.ImageWidth;% check which points in the subimage 2 are contained in the initial image 2
675    check2_y=subrange2_y>=1 & subrange2_y<=par_civ.ImageHeight;
676   
677    image1_crop(check1_y,check1_x)=par_civ.ImageA(subrange1_y(check1_y),subrange1_x(check1_x));%extract a subimage (correlation box) from image A
678    image2_crop(check2_y,check2_x)=par_civ.ImageB(subrange2_y(check2_y),subrange2_x(check2_x));%extract a larger subimage (search box) from image B
679    image1_mean=mean(mean(image1_crop));
680    image2_mean=mean(mean(image2_crop));
681    %threshold on image minimum
682    if check_MinIma && (image1_mean < par_civ.MinIma || image2_mean < par_civ.MinIma)
683        F(ivec)=3;
684    end
685    %threshold on image maximum
686    if check_MaxIma && (image1_mean > par_civ.MaxIma || image2_mean > par_civ.MaxIma)
687        F(ivec)=3;
688    end
689    %         end
690    if F(ivec)~=3
691        image1_crop=image1_crop-image1_mean;%substract the mean
692        image2_crop=image2_crop-image2_mean;
693        if CheckDecimal
694            xi=(1:mesh:size(image1_crop,2));
695            yi=(1:mesh:size(image1_crop,1))';
696            if CheckDeformation
697                [XI,YI]=meshgrid(xi-ceil(size(image1_crop,2)/2),yi-ceil(size(image1_crop,1)/2));
698                XIant=XI-par_civ.DUDX(ivec)*XI-par_civ.DUDY(ivec)*YI+ceil(size(image1_crop,2)/2);
699                YIant=YI-par_civ.DVDX(ivec)*XI-par_civ.DVDY(ivec)*YI+ceil(size(image1_crop,1)/2);
700                image1_crop=interp2(image1_crop,XIant,YIant);
701            else
702                image1_crop=interp2(image1_crop,xi,yi);
703            end
704            xi=(1:mesh:size(image2_crop,2));
705            yi=(1:mesh:size(image2_crop,1))';
706            image2_crop=interp2(image2_crop,xi,yi);
707        end
708        sum_square=sum(sum(image1_crop.*image1_crop));
709        %reference: Oliver Pust, PIV: Direct Cross-Correlation
710        result_conv= conv2(image2_crop,flipdim(flipdim(image1_crop,2),1),'valid');
711        corrmax= max(max(result_conv));
712        result_conv=(result_conv/corrmax)*255; %normalize, peak=always 255
713        %Find the correlation max, at 255
714        [y,x] = find(result_conv==255,1);
715        if ~isempty(y) && ~isempty(x)
716            try
717                if par_civ.CorrSmooth==1
718                    [vector,F(ivec)] = SUBPIXGAUSS (result_conv,x,y);
719                elseif par_civ.CorrSmooth==2
720                    [vector,F(ivec)] = SUBPIX2DGAUSS (result_conv,x,y);
721                end
722                utable(ivec)=vector(1)*mesh+shiftx(ivec);
723                vtable(ivec)=vector(2)*mesh+shifty(ivec);
724                xtable(ivec)=iref+utable(ivec)/2-0.5;% convec flow (velocity taken at the point middle from imgae 1 and 2)
725                ytable(ivec)=jref+vtable(ivec)/2-0.5;% and position of pixel 1=0.5 (convention for image coordinates=0 at the edge)
726                iref=round(xtable(ivec));% image index for the middle of the vector
727                jref=round(ytable(ivec));
728                if checkmask && par_civ.Mask(jref,iref)<200 && par_civ.Mask(jref,iref)>=100
729                    utable(ivec)=0;
730                    vtable(ivec)=0;
731                    F(ivec)=3;
732                end
733                ctable(ivec)=corrmax/sum_square;% correlation value
734            catch ME
735                F(ivec)=3;
736            end
737        else
738            F(ivec)=3;
739        end
740    end
741end
742result_conv=result_conv*corrmax/(255*sum_square);% keep the last correlation matrix for output
743
744%------------------------------------------------------------------------
745% --- Find the maximum of the correlation function after interpolation
746function [vector,F] = SUBPIXGAUSS (result_conv,x,y)
747%------------------------------------------------------------------------
748vector=[0 0]; %default
749F=0;
750[npy,npx]=size(result_conv);
751
752% if (x <= (size(result_conv,1)-1)) && (y <= (size(result_conv,1)-1)) && (x >= 1) && (y >= 1)
753    %the following 8 lines are copyright (c) 1998, Uri Shavit, Roi Gurka, Alex Liberzon, Technion ï¿œ Israel Institute of Technology
754    %http://urapiv.wordpress.com
755    peaky = y;
756    if y <= npy-1 && y >= 1
757        f0 = log(result_conv(y,x));
758        f1 = real(log(result_conv(y-1,x)));
759        f2 = real(log(result_conv(y+1,x)));
760        peaky = peaky+ (f1-f2)/(2*f1-4*f0+2*f2);
761    else
762        F=-2; % warning flag for vector truncated by the limited search box
763    end
764    peakx=x;
765    if x <= npx-1 && x >= 1
766        f0 = log(result_conv(y,x));
767        f1 = real(log(result_conv(y,x-1)));
768        f2 = real(log(result_conv(y,x+1)));
769        peakx = peakx+ (f1-f2)/(2*f1-4*f0+2*f2);
770    else
771        F=-2; % warning flag for vector truncated by the limited search box
772    end
773    vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
774
775%------------------------------------------------------------------------
776% --- Find the maximum of the correlation function after interpolation
777function [vector,F] = SUBPIX2DGAUSS (result_conv,x,y)
778%------------------------------------------------------------------------
779vector=[0 0]; %default
780F=-2;
781peaky=y;
782peakx=x;
783[npy,npx]=size(result_conv);
784if (x <= npx-1) && (y <= npy-1) && (x >= 1) && (y >= 1)
785    F=0;
786    for i=-1:1
787        for j=-1:1
788            %following 15 lines based on
789            %H. Nobach ï¿œ M. Honkanen (2005)
790            %Two-dimensional Gaussian regression for sub-pixel displacement
791            %estimation in particle image velocimetry or particle position
792            %estimation in particle tracking velocimetry
793            %Experiments in Fluids (2005) 38: 511ï¿œ515
794            c10(j+2,i+2)=i*log(result_conv(y+j, x+i));
795            c01(j+2,i+2)=j*log(result_conv(y+j, x+i));
796            c11(j+2,i+2)=i*j*log(result_conv(y+j, x+i));
797            c20(j+2,i+2)=(3*i^2-2)*log(result_conv(y+j, x+i));
798            c02(j+2,i+2)=(3*j^2-2)*log(result_conv(y+j, x+i));
799        end
800    end
801    c10=(1/6)*sum(sum(c10));
802    c01=(1/6)*sum(sum(c01));
803    c11=(1/4)*sum(sum(c11));
804    c20=(1/6)*sum(sum(c20));
805    c02=(1/6)*sum(sum(c02));
806    deltax=(c11*c01-2*c10*c02)/(4*c20*c02-c11^2);
807    deltay=(c11*c10-2*c01*c20)/(4*c20*c02-c11^2);
808    if abs(deltax)<1
809        peakx=x+deltax;
810    end
811    if abs(deltay)<1
812        peaky=y+deltay;
813    end
814end
815vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
816
817%'RUN_FIX': function for fixing velocity fields:
818%-----------------------------------------------
819% RUN_FIX(filename,field,flagindex,thresh_vecC,thresh_vel,iter,flag_mask,maskname,fileref,fieldref)
820%
821%filename: name of the netcdf file (used as input and output)
822%field: structure specifying the names of the fields to fix (depending on civ1 or civ2)
823    %.vel_type='civ1' or 'civ2';
824    %.nb=name of the dimension common to the field to fix ('nb_vectors' for civ1);
825    %.fixflag=name of fix flag variable ('vec_FixFlag' for civ1)
826%flagindex: flag specifying which values of vec_f are removed:
827        % if flagindex(1)=1: vec_f=-2 vectors are removed
828        % if flagindex(2)=1: vec_f=3 vectors are removed
829        % if flagindex(3)=1: vec_f=2 vectors are removed (if iter=1) or vec_f=4 vectors are removed (if iter=2)
830%iter=1 for civ1 fields and iter=2 for civ2 fields
831%thresh_vecC: threshold in the image correlation vec_C
832%flag_mask: =1 mask used to remove vectors (0 else)
833%maskname: name of the mask image file for fix
834%thresh_vel: threshold on velocity, or on the difference with the reference file fileref if exists
835%inf_sup=1: remove values smaller than threshold thresh_vel, =2, larger than threshold
836%fileref: .nc file name for a reference velocity (='': refrence 0 used)
837%fieldref: 'civ1','filter1'...feld used in fileref
838
839function FF=fix(Param,F,C,U,V,X,Y)
840FF=zeros(size(F));%default
841
842%criterium on warn flags
843FlagName={'CheckFmin2','CheckF2','CheckF3','CheckF4'};
844FlagVal=[-2 2 3 4];
845for iflag=1:numel(FlagName)
846    if isfield(Param,FlagName{iflag}) && Param.(FlagName{iflag})
847        FF=(FF==1| F==FlagVal(iflag));
848    end
849end
850%criterium on correlation values
851if isfield (Param,'MinCorr')
852    FF=FF==1 | C<Param.MinCorr;
853end
854if (isfield(Param,'MinVel')&&~isempty(Param.MinVel))||(isfield (Param,'MaxVel')&&~isempty(Param.MaxVel))
855    Umod= U.*U+V.*V;
856    if isfield (Param,'MinVel')&&~isempty(Param.MinVel)
857        FF=FF==1 | Umod<(Param.MinVel*Param.MinVel);
858    end
859    if isfield (Param,'MaxVel')&&~isempty(Param.MaxVel)
860        FF=FF==1 | Umod>(Param.MaxVel*Param.MaxVel);
861    end
862end
863
864
865%------------------------------------------------------------------------
866% --- determine the list of index pairs of processing file
867function [i1_series,i2_series,j1_series,j2_series,check_bounds,NomTypeNc]=...
868    find_pair_indices(str_civ,i_series,j_series,MinIndex,MaxIndex)
869%------------------------------------------------------------------------
870i1_series=i_series;% set of first image numbers
871i2_series=i_series;
872j1_series=ones(size(i_series));% set of first image numbers
873j2_series=ones(size(i_series));
874check_bounds=false(size(i_series));
875r=regexp(str_civ,'^\D(?<ind>[i|j])= (?<num1>\d+)\|(?<num2>\d+)','names');
876if ~isempty(r)
877    mode=['D' r.ind];
878    ind1=str2num(r.num1);
879    ind2=str2num(r.num2);
880else
881    mode='burst';
882    r=regexp(str_civ,'^j= (?<num1>[a-z])-(?<num2>[a-z])','names');
883    if ~isempty(r)
884        NomTypeNc='_1ab';
885    else
886        r=regexp(str_civ,'^j= (?<num1>[A-Z])-(?<num2>[A-Z])','names');
887        if ~isempty(r)
888            NomTypeNc='_1AB';
889        else
890            r=regexp(str_civ,'^j= (?<num1>\d+)-(?<num2>\d+)','names');
891            if ~isempty(r)
892                NomTypeNc='_1_1-2';
893            end           
894        end
895    end
896    if isempty(r)
897        display('wrong pair mode input option')
898    else
899    ind1=stra2num(r.num1);
900    ind2=stra2num(r.num2);
901    end
902end
903if strcmp (mode,'Di')
904    i1_series=i_series-ind1;% set of first image numbers
905    i2_series=i_series+ind2;
906    check_bounds=i1_series<MinIndex(1,1) | i2_series>MaxIndex(1,1);
907    if isempty(j_series)
908        NomTypeNc='_1-2';
909    else
910        j1_series=j_series;
911        j2_series=j_series;
912        NomTypeNc='_1-2_1';
913    end
914elseif strcmp (mode,'Dj')
915    j1_series=j_series-ind1;
916    j2_series=j_series+ind2;
917    check_bounds=j1_series<MinIndex(1,2) | j2_series>MaxIndex(1,2);
918    NomTypeNc='_1_1-2';
919else  %bursts
920    j1_series=ind1*ones(size(i_series));
921    j2_series=ind2*ones(size(i_series));
922end
923
924%     if length(indsel)>=1
925%         firstind=indsel(1);
926%         lastind=indsel(end);
927%         set(handles.first_j,'String',num2str(ref_j(firstind)))%update the display of first and last fields
928%         set(handles.last_j,'String',num2str(ref_j(lastind)))
929%         ref_j=ref_j(indsel);
930%         j1_civ1=j1_civ1(indsel);
931%         j2_civ1=j2_civ1(indsel);
932%         j1_civ2=j1_civ2(indsel);
933%         j2_civ2=j2_civ2(indsel);
934%     end
935% elseif isequal(mode,'pair j1-j2') %case of bursts (png_old or png_2D)
936%     displ_num=get(handles.ListPairCiv1,'UserData');
937%     i1_civ1=ref_i;
938%     i2_civ1=ref_i;
939%     j1_civ1=displ_num(1,index_civ1);
940%     j2_civ1=displ_num(2,index_civ1);
941%     i1_civ2=ref_i;
942%     i2_civ2=ref_i;
943%     j1_civ2=displ_num(1,index_civ2);
944%     j2_civ2=displ_num(2,index_civ2);
945% elseif isequal(mode,'displacement')
946%     i1_civ1=ref_i;
947%     i2_civ1=ref_i;
948%     j1_civ1=ref_j;
949%     j2_civ1=ref_j;
950%     i1_civ2=ref_i;
951%     i2_civ2=ref_i;
952%     j1_civ2=ref_j;
953%     j2_civ2=ref_j;
954% end
955
956
957
958
Note: See TracBrowser for help on using the repository browser.