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

Last change on this file since 1156 was 1156, checked in by sommeria, 4 months ago

bugs repaired

File size: 65.4 KB
Line 
1%'civ_series': PIV function activated by the general GUI series
2% --- call the sub-functions:
3%   civ: PIV function itself
4%   detect_false: put a flag to false vectors after detection by various criteria
5%   filter_tps: make interpolation-smoothing
6%------------------------------------------------------------------------
7% function [Data,errormsg,result_conv]= civ_series(Param)
8%
9%OUTPUT
10% Data=structure containing the PIV results and information on the processing parameters
11% errormsg=error message char string, decd ..fault=''
12% resul_conv: image inter-correlation function for the last grid point (used for tests)
13%
14%INPUT:
15% Param: Matlab structure of input  parameters
16%     Param contains info of the GUI series using the fct read_GUI.
17%     Param.Action.RUN = 0 (to set the status of the GUI series) or =1 to RUN the computation
18%     Param.InputTable: sets the input file(s)
19%           if absent, the fct looks for input data in Param.ActionInput     (test mode)
20%     Param.OutputSubDir: sets the folder name of output file(s,
21%           if absent no file is produced, result in the output structure Data (test mode)
22%     Param.ActionInput: substructure with the parameters provided by the GUI civ_input
23%                      .Civ1: parameters for civ1cc
24%                      .Fix1: parameters for detect_false1
25%                      .Patch1:
26%                      .Civ2: for civ2
27%                      .Fix2:
28%                      .Patch2:
29
30%=======================================================================
31% Copyright 2008-2024, LEGI UMR 5519 / CNRS UGA G-INP, Grenoble, France
32%   http://www.legi.grenoble-inp.fr
33%   Joel.Sommeria - Joel.Sommeria (A) univ-grenoble-alpes.fr
34%
35%     This file is part of the toolbox UVMAT.
36%
37%     UVMAT is free software; you can redistribute it and/or modify
38%     it under the terms of the GNU General Public License as published
39%     by the Free Software Foundation; either version 2 of the license,
40%     or (at your option) any later version.
41%
42%     UVMAT is distributed in the hope that it will be useful,
43%     but WITHOUT ANY WARRANTY; without even the implied warranty of
44%     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
45%     GNU General Public License (see LICENSE.txt) for more details.
46%=======================================================================
47
48function [Data,errormsg,result_conv]= civ_series(Param)
49errormsg='';
50
51%% set the input elements needed on the GUI series when the action is selected in the menu ActionName or InputTable refreshed
52if isstruct(Param) && isequal(Param.Action.RUN,0)% function activated from the GUI series but not RUN
53    path_series=fileparts(which('series'));
54    addpath(fullfile(path_series,'series'))
55    Data=civ_input(Param);% introduce the civ parameters using the GUI civ_input
56    % TODO: change from guide to App: modify the input procedure, adapt read_GUI function
57    %App=civ_input_App
58    %Data=civ_input_App(Param);% introduce the civ parameters using the GUI civ_input
59    % if isempty(App)
60    %     Data=Param;% if  civ_input has been cancelled, keep previous parameters
61    % end
62    Data.Program=mfilename;%gives the name of the current function
63    Data.AllowInputSort='off';% allow alphabetic sorting of the list of input file SubDir (options 'off'/'on', 'off' by default)
64    Data.WholeIndexRange='off';% prescribes the file index ranges from min to max (options 'off'/'on', 'off' by default)
65    Data.NbSlice='off'; %nbre of slices ('off' by default)
66    Data.VelType='off';% menu for selecting the velocity type (options 'off'/'one'/'two',  'off' by default)
67    Data.FieldName='on';% menu for selecting the field (s) in the input file(options 'off'/'one'/'two', 'off' by default)
68    Data.FieldTransform = 'off';%can use a transform function
69    Data.ProjObject='off';%can use projection object(option 'off'/'on',
70    Data.Mask='off';%can use mask option   (option 'off'/'on', 'off' by default)
71    Data.OutputDirExt='.civ';%set the output dir extension
72    Data.OutputSubDirMode='last'; %select the last subDir in the input table as root of the output subdir name (option 'all'/'first'/'last', 'all' by default)
73    Data.OutputFileMode='NbInput_i';% one output file expected per value of i index (used for waitbar)
74    Data.CheckOverwriteVisible='on'; % manage the overwrite of existing files (default=1)
75    if isfield(Data,'ActionInput') && isfield(Data.ActionInput,'PairIndices') && strcmp(Data.ActionInput.PairIndices.ListPairMode,'pair j1-j2')
76        Data.IndexRange_j='off';%no j index display in series
77        % if isfield(Data.ActionInput.PairIndices,'ListPairCiv2')
78        %     str_civ=Data.ActionInput.PairIndices.ListPairCiv2;
79        % else
80        %     str_civ=Data.ActionInput.PairIndices.ListPairCiv1;
81        % end
82        % r=regexp(str_civ,'^j= (?<num1>[a-z])-(?<num2>[a-z])','names');
83        % if isempty(r)
84        %     r=regexp(str_civ,'^j= (?<num1>[A-Z])-(?<num2>[A-Z])','names');
85        %     if isempty(r)
86        %         r=regexp(str_civ,'^j= (?<num1>\d+)-(?<num2>\d+)','names');
87        %     end
88        % end
89        % if ~isempty(r)
90        %     Data.j_index_1=stra2num(r.num1);
91        %     Data.j_index_2=stra2num(r.num2);
92        % end
93    else
94        Data.IndexRange_j='on';% j index display in series if relevant
95    end
96    return
97end
98
99%% read input parameters from an xml file if input is a file name (batch mode)
100checkrun=1;
101if ischar(Param)
102    Param=xml2struct(Param);% read Param as input file (batch case)
103    checkrun=0;
104end
105
106%% test input
107if ~isfield(Param,'ActionInput')
108    disp_uvmat('ERROR','no parameter set for PIV',checkrun)
109    return
110end
111iview_A=0;%default values
112NbField=1;
113RUNHandle=[];
114CheckInputFile=isfield(Param,'InputTable');%= 1 in test use for TestCiv (no nc file involved)
115CheckOutputFile=isfield(Param,'OutputSubDir');%= 1 in test use for TestPatch (no nc file produced)
116
117%% input files and indexing (skipped in Test mode)
118if CheckInputFile
119    hseries=findobj(allchild(0),'Tag','series');
120    RUNHandle=findobj(hseries,'Tag','RUN');%handle of RUN button in GUI series
121    WaitbarHandle=findobj(hseries,'Tag','Waitbar');%handle of waitbar in GUI series
122    MaxIndex_i=Param.IndexRange.MaxIndex_i;
123    MinIndex_i=Param.IndexRange.MinIndex_i;
124    MaxIndex_j=ones(size(MaxIndex_i));MinIndex_j=ones(size(MinIndex_i));
125    if isfield(Param.IndexRange,'MaxIndex_j')&& isfield(Param.IndexRange,'MinIndex_j')
126        MaxIndex_j=Param.IndexRange.MaxIndex_j;
127        MinIndex_j=Param.IndexRange.MinIndex_j;
128    end
129    if isfield(Param,'InputTable')
130        [~,i1_series,i2_series,j1_series,j2_series]=get_file_series(Param);
131        iview_A=0;% series index (iview) for the first image series
132        iview_B=0;% series index (iview) for the second image series (only non zero for option 'shift' comparing two image series )
133        if Param.ActionInput.CheckCiv1
134            iview_A=1;% usual PIV, the image series is on the first line of the table
135        elseif Param.ActionInput.CheckCiv2 % civ2 is performed without Civ1, a netcdf file series is needed in the first table line
136            iview_A=2;% the second line is used for the input images of Civ2
137        end
138        if iview_A~=0
139            RootPath_A=Param.InputTable{iview_A,1};
140            RootFile_A=Param.InputTable{iview_A,3};
141            SubDir_A=Param.InputTable{iview_A,2};
142            NomType_A=Param.InputTable{iview_A,4};
143            FileExt_A=Param.InputTable{iview_A,5};
144            if iview_B==0
145                iview_B=iview_A;% the second image series is the same as the first
146            end
147            RootPath_B=Param.InputTable{iview_B,1};
148            RootFile_B=Param.InputTable{iview_B,3};
149            SubDir_B=Param.InputTable{iview_B,2};
150            NomType_B=Param.InputTable{iview_B,4};
151            FileExt_B=Param.InputTable{iview_B,5};
152        end
153       
154        PairCiv2='';
155        switch Param.ActionInput.ListCompareMode
156            case 'PIV'
157                PairCiv1=Param.ActionInput.PairIndices.ListPairCiv1;
158                if isfield(Param.ActionInput.PairIndices,'ListPairCiv2')
159                    PairCiv2=Param.ActionInput.PairIndices.ListPairCiv2;%string which determines the civ2 pair
160                end
161                if iview_A==1% if Civ1 is performed
162                    [i1_series_Civ1,i2_series_Civ1,j1_series_Civ1,j2_series_Civ1,check_bounds,NomTypeNc]=...
163                        find_pair_indices(PairCiv1,i1_series{1},j1_series{1},MinIndex_i,MaxIndex_i,MinIndex_j,MaxIndex_j);
164                    if ~isempty(PairCiv2)
165                        [i1_series_Civ2,i2_series_Civ2,j1_series_Civ2,j2_series_Civ2,check_bounds_Civ2]=...
166                            find_pair_indices(PairCiv2,i1_series{1},j1_series{1},MinIndex_i(1),MaxIndex_i(1),MinIndex_j(1),MaxIndex_j(1));
167                        check_bounds=check_bounds | check_bounds_Civ2;
168                    end
169                else% we start from an existing Civ1 file
170                    i1_series_Civ1=i1_series{1};
171                    i2_series_Civ1=i2_series{1};
172                    j1_series_Civ1=j1_series{1};
173                    j2_series_Civ1=j2_series{1};
174                    NomTypeNc=Param.InputTable{1,4};
175                    if ~isempty(PairCiv2)
176                        [i1_series_Civ2,i2_series_Civ2,j1_series_Civ2,j2_series_Civ2,check_bounds,NomTypeNc]=...
177                            find_pair_indices(PairCiv2,i1_series{2},j1_series{2},MinIndex_i(2),MaxIndex_i(2),MinIndex_j(2),MaxIndex_j(2));
178                    end
179                end
180            case 'displacement'
181                if isfield(Param.ActionInput,'OriginIndex')
182                i1_series_Civ1=Param.ActionInput.OriginIndex*ones(size(i1_series{1}));
183                else
184                    i1_series_Civ1=ones(size(i1_series{1}));
185                end
186                i1_series_Civ2=i1_series_Civ1;
187                i2_series_Civ1=i1_series{1};
188                i2_series_Civ2=i1_series{1};
189                j1_series_Civ1=[];% no j index variation for the ref image
190                j1_series_Civ2=[];
191                if isempty(j1_series{1})
192                    j2_series_Civ1=ones(size(i1_series_Civ1));
193                else
194                    j2_series_Civ1=j1_series{1};% if j index exist 
195                end
196                j2_series_Civ2=j2_series_Civ1;
197                NomTypeNc='_1';
198        end
199        %determine frame indices for input with movie or other multiframe input file
200        if isempty(j1_series_Civ1)% simple movie with index i
201            FrameIndex_A_Civ1=i1_series_Civ1;
202            FrameIndex_B_Civ1=i2_series_Civ1;
203            j1_series_Civ1=ones(size(i1_series_Civ1));
204            if strcmp(Param.ActionInput.ListCompareMode,'PIV')
205            j2_series_Civ1=ones(size(i1_series_Civ1));
206            end
207        else % movie for each burst or volume (index j)
208            FrameIndex_A_Civ1=j1_series_Civ1;
209            FrameIndex_B_Civ1=j2_series_Civ1;
210        end
211        if isempty(PairCiv2)
212            FrameIndex_A_Civ2=FrameIndex_A_Civ1;
213            FrameIndex_B_Civ2=FrameIndex_B_Civ1;
214        else
215            if isempty(j1_series_Civ2)
216                FrameIndex_A_Civ2=i1_series_Civ2;
217                FrameIndex_B_Civ2=i2_series_Civ2;
218                j1_series_Civ2=ones(size(i1_series_Civ2));
219                if strcmp(Param.ActionInput.ListCompareMode,'PIV')
220                j2_series_Civ2=ones(size(i1_series_Civ2));
221                end
222            else
223                FrameIndex_A_Civ2=j1_series_Civ2;
224                FrameIndex_B_Civ2=j2_series_Civ2;
225            end
226        end
227        if isempty(i1_series_Civ1)||(~isempty(PairCiv2) && isempty(i1_series_Civ2))
228            disp_uvmat('ERROR','no image pair for civ in the input file index range',checkrun)
229            return
230        end
231    end
232   
233    %% check the first image pair
234        if Param.ActionInput.CheckCiv1% Civ1 is performed
235            NbField=numel(i1_series_Civ1);
236        elseif Param.ActionInput.CheckCiv2 % Civ2 is performed without Civ1
237            NbField=numel(i1_series_Civ2);
238        else
239            NbField=numel(i1_series_Civ1);% no image used (only detect_false or patch) TO CHECK
240        end
241
242    %% Output directory
243    OutputDir='';
244    if CheckOutputFile
245        OutputDir=[Param.OutputSubDir Param.OutputDirExt];
246    end
247end
248
249%% prepare output Data
250ListGlobalAttribute={'Conventions','Program','CivStage'};
251Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
252Data.Program='civ_series';
253if isfield(Param,'UvmatRevision')
254    Data.Program=[Data.Program ', uvmat r' Param.UvmatRevision];
255end
256Data.CivStage=0;%default
257
258%% get timing from the ImaDoc file or input video
259if iview_A~=0
260    XmlFileName=find_imadoc(RootPath_A,SubDir_A);
261    Time=[];
262    if ~isempty(XmlFileName)
263        XmlData=imadoc2struct(XmlFileName);
264        if isfield(XmlData,'Time')
265            Time=XmlData.Time;
266        end
267        if isfield(XmlData,'Camera')
268            if isfield(XmlData.Camera,'NbSlice')&& ~isempty(XmlData.Camera.NbSlice)
269                NbSlice_calib{iview}=XmlData.Camera.NbSlice;% Nbre of slices for Zindex in phys transform
270                if ~isequal(NbSlice_calib{iview},NbSlice_calib{1})
271                    msgbox_uvmat('WARNING','inconsistent number of Z indices for the two field series');
272                end
273            end
274            if isfield(XmlData.Camera,'TimeUnit')&& ~isempty(XmlData.Camera.TimeUnit)
275                TimeUnit=XmlData.Camera.TimeUnit;
276            end
277        end
278    end
279end
280
281%%%%% MAIN LOOP %%%%%%
282maskoldname='';% initiate the mask name
283FileType_A='';
284FileType_B='';
285CheckOverwrite=1;%default
286if isfield(Param,'CheckOverwrite')
287    CheckOverwrite=Param.CheckOverwrite;
288end
289for ifield=1:NbField
290    tstart=tic;
291    time_civ1=0;
292  time_patch1=0;
293  time_civ2=0;
294  time_patch2=0;
295    if ~isempty(RUNHandle)% update the waitbar in interactive mode with GUI series  (checkrun=1)
296        update_waitbar(WaitbarHandle,ifield/NbField)
297        if  checkrun && ~strcmp(get(RUNHandle,'BusyAction'),'queue')
298            disp('program stopped by user')
299            break
300        end
301    end
302    if CheckInputFile
303        OutputPath=fullfile(Param.OutputPath,num2str(Param.Experiment),num2str(Param.Device));
304        if iview_A==0 % no nc file has been entered
305            ncfile=fullfile_uvmat(OutputPath,Param.InputTable{1,2},Param.InputTable{1,3},Param.InputTable{1,5},...
306                NomTypeNc,i1_series_Civ1(ifield),i2_series_Civ1(ifield),j1_series_Civ1(ifield),j2_series_Civ1(ifield));
307        else% an existing nc file has been entered
308            if iview_A==1% if Civ1 is performed
309                Civ1Dir=OutputDir;
310            else
311                Civ1Dir=Param.InputTable{1,2};
312            end
313            if strcmp(Param.ActionInput.ListCompareMode,'PIV')
314                ncfile=fullfile_uvmat(OutputPath,Civ1Dir,RootFile_A,'.nc',NomTypeNc,i1_series_Civ1(ifield),i2_series_Civ1(ifield),...
315                    j1_series_Civ1(ifield),j2_series_Civ1(ifield));
316            else
317                ncfile=fullfile_uvmat(OutputPath,Civ1Dir,RootFile_A,'.nc',NomTypeNc,i2_series_Civ1(ifield),[],...
318                    j1_series_Civ1(ifield),j2_series_Civ1(ifield));
319            end
320        end
321        ncfile_out=ncfile;% by default
322       
323        if isfield (Param.ActionInput,'Civ2')
324            i1_civ2=i1_series_Civ2(ifield);
325            i2_civ2=i1_civ2;
326            if ~isempty(i2_series_Civ2)
327                i2_civ2=i2_series_Civ2(ifield);
328            end
329            j1_civ2=1;
330            if ~isempty(j1_series_Civ2)
331                j1_civ2=j1_series_Civ2(ifield);
332            end
333            j2_civ2=i1_civ2;
334            if ~isempty(j2_series_Civ2)
335                j2_civ2=j2_series_Civ2(ifield);
336            end
337            if strcmp(Param.ActionInput.ListCompareMode,'PIV')
338                ncfile_out=fullfile_uvmat(OutputPath,OutputDir,RootFile_A,'.nc',NomTypeNc,i1_civ2,i2_civ2,j1_civ2,j2_civ2);
339            else % displacement
340                ncfile_out=fullfile_uvmat(OutputPath,OutputDir,RootFile_A,'.nc',NomTypeNc,i2_civ2,[],j2_civ2);
341            end
342        end
343        if ~CheckOverwrite && exist(ncfile_out,'file')
344            disp(['existing output file ' ncfile_out ' already exists, skip to next field'])
345            continue% skip iteration if the mode overwrite is desactivated and the result file already exists
346        end
347    end
348    ImageName_A='';ImageName_B='';%default
349    VideoObject_A=[];VideoObject_B=[];
350   
351    %% Civ1
352    % if Civ1 computation is requested
353    if isfield (Param.ActionInput,'Civ1')
354        if CheckInputFile
355            disp('civ1 started')
356        end
357        par_civ1=Param.ActionInput.Civ1;% parameters for civ1
358        if CheckInputFile % read input images (except in mode Test where it is introduced directly in Param.ActionInput.Civ1.ImageNameA and B)
359            try
360                if strcmp(Param.ActionInput.ListCompareMode,'displacement')
361                    ImageName_A=Param.ActionInput.RefFile;
362                else
363                ImageName_A=fullfile_uvmat(RootPath_A,SubDir_A,RootFile_A,FileExt_A,NomType_A,i1_series_Civ1(ifield),[],j1_series_Civ1(ifield));
364                end
365                if strcmp(FileExt_A,'.nc')% case of input images in format netcdf
366                    FieldName_A=Param.InputFields.FieldName;
367                    [DataIn,~,~,errormsg]=nc2struct(ImageName_A,{FieldName_A});
368                    par_civ1.ImageA=DataIn.(FieldName_A);
369                else % usual image formats for image A
370                    if isempty(FileType_A)% open the image object if not already done in case of movie input
371                        [FileInfo_A,VideoObject_A]=get_file_info(ImageName_A);
372                        FileType_A=FileInfo_A.FileType;
373                        if isempty(Time) && ~isempty(find(strcmp(FileType_A,{'mmreader','video','cine_phantom'}), 1))% case of video input
374                            Time=zeros(FileInfo_A.NumberOfFrames+1,2);
375                            Time(:,2)=(0:1/FileInfo_A.FrameRate:(FileInfo_A.NumberOfFrames)/FileInfo_A.FrameRate)';
376                        end
377                        if ~isempty(FileType_A) && isempty(Time)% Time = index i +0.001 index j by default
378                            MaxIndex_i=max(i2_series_Civ1);
379                            MaxIndex_j=max(j2_series_Civ1);
380                            Time=(1:MaxIndex_i)'*ones(1,MaxIndex_j);
381                            Time=Time+0.001*ones(MaxIndex_i,1)*(1:MaxIndex_j);
382                            Time=[zeros(1,MaxIndex_j);Time];% insert a first line of zeros
383                            Time=[zeros(MaxIndex_i+1,1) Time];% insert a first column of zeros
384                        end
385                    end
386                    if isempty(regexp(ImageName_A,'(^http://)|(^https://)', 'once')) && ~exist(ImageName_A,'file')
387                        disp([ImageName_A ' missing'])
388                        continue
389                    end
390                    tsart_input=tic;
391                    [par_civ1.ImageA,VideoObject_A] = read_image(ImageName_A,FileType_A,VideoObject_A,FrameIndex_A_Civ1(ifield));
392                    time_input=toc(tsart_input);
393                end
394                ImageName_B=fullfile_uvmat(RootPath_B,SubDir_B,RootFile_B,FileExt_B,NomType_B,i2_series_Civ1(ifield),[],j2_series_Civ1(ifield));
395                if strcmp(FileExt_B,'.nc') % case of input images in format netcdf
396                    FieldName_B=Param.InputFields.FieldName;
397                    [DataIn,~,~,errormsg]=nc2struct(ImageName_B,{FieldName_B});
398                    par_civ1.ImageB=DataIn.(FieldName_B);
399                else % usual image formats for image B
400                    if isempty(FileType_B)
401                        [FileInfo_B,VideoObject_B]=get_file_info(ImageName_B);
402                        FileType_B=FileInfo_B.FileType;
403                    end
404                    if isempty(regexp(ImageName_B,'(^http://)|(^https://)', 'once')) && ~exist(ImageName_B,'file')
405                        disp([ImageName_B ' missing'])
406                        continue
407                    end
408                    [par_civ1.ImageB,VideoObject_B] = read_image(ImageName_B,FileType_B,VideoObject_B,FrameIndex_B_Civ1(ifield));
409                end
410            catch ME % display errors in reading input images
411                if ~isempty(ME.message)
412                    disp_uvmat('ERROR', ['error reading input image: ' ME.message],checkrun)
413                    continue
414                end
415            end
416            par_civ1.ImageWidth=size(par_civ1.ImageA,2);
417            par_civ1.ImageHeight=size(par_civ1.ImageA,1);
418            list_param=(fieldnames(Param.ActionInput.Civ1))';
419            list_param(strcmp('TestCiv1',list_param))=[];% remove the parameter TestCiv1 from the list
420            Civ1_param=regexprep(list_param,'^.+','Civ1_$0');% insert 'Civ1_' before  each string in list_param
421            Civ1_param=[{'Civ1_ImageA','Civ1_ImageB','Civ1_Time','Civ1_Dt'} Civ1_param]; %insert the names of the two input images
422            %indicate the values of all the global attributes in the output data
423            Data.Civ1_ImageA=ImageName_A;
424            Data.Civ1_ImageB=ImageName_B;
425            i1=i1_series_Civ1(ifield);
426            i2=i1;
427            if ~isempty(i2_series_Civ1)
428                i2=i2_series_Civ1(ifield);
429            end
430            j1=1;
431            if ~isempty(j1_series_Civ1)
432                j1=j1_series_Civ1(ifield);
433            end
434            j2=j1;
435            if ~isempty(j2_series_Civ1)
436                j2=j2_series_Civ1(ifield);
437            end
438            if strcmp(Param.ActionInput.ListCompareMode,'displacement')
439                Data.Civ1_Time=Time(i2+1,j2+1);% the Time is the Time of the secodn image
440                Data.Civ1_Dt=1;% Time interval is 1, to yield displacement instead of velocity=displacement/Dt at reading
441            else
442                Data.Civ1_Time=(Time(i2+1,j2+1)+Time(i1+1,j1+1))/2;% the Time is the Time at the middle of the image pair
443                Data.Civ1_Dt=Time(i2+1,j2+1)-Time(i1+1,j1+1);
444            end
445            for ilist=1:length(list_param)
446                Data.(Civ1_param{4+ilist})=Param.ActionInput.Civ1.(list_param{ilist});
447            end
448            Data.ListGlobalAttribute=[ListGlobalAttribute Civ1_param];
449            Data.CivStage=1;
450        else
451            i1=Param.ActionInput.PairIndices.ref_i; %case of TESTmode
452        end
453        % set the list of variables
454        Data.ListVarName={'Civ1_X','Civ1_Y','Civ1_U','Civ1_V','Civ1_C','Civ1_FF'};%  cell array containing the names of the fields to record
455        Data.VarDimName={'nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1','nb_vec_1'};
456        Data.VarAttribute{1}.Role='coord_x';
457        Data.VarAttribute{2}.Role='coord_y';
458        Data.VarAttribute{3}.Role='vector_x';
459        Data.VarAttribute{4}.Role='vector_y';
460        Data.VarAttribute{5}.Role='ancillary';
461        Data.VarAttribute{6}.Role='errorflag';
462        % case of mask
463        if par_civ1.CheckMask&&~isempty(par_civ1.Mask)
464            if isfield(par_civ1,'NbSlice')
465                [RootPath_mask,SubDir_mask,RootFile_mask,~,~,~,~,Ext_mask]=fileparts_uvmat(Param.ActionInput.Civ1.Mask);
466                i1_mask=mod(i1-1,par_civ1.NbSlice)+1;
467                maskname=fullfile_uvmat(RootPath_mask,SubDir_mask,RootFile_mask,Ext_mask,'_1',i1_mask);
468            else
469                maskname=Param.ActionInput.Civ1.Mask;
470            end
471            if strcmp(maskoldname,maskname)% mask exist, not already read in civ1
472                par_civ1.Mask=mask; %use mask already opened
473            else
474                if ~isempty(regexp(maskname,'(^http://)|(^https://)', 'once'))|| exist(maskname,'file')
475                    try
476                        par_civ1.Mask=imread(maskname);%update the mask, an store it for future use
477                    catch ME
478                        if ~isempty(ME.message)
479                            errormsg=['error reading input image: ' ME.message];
480                            disp_uvmat('ERROR',errormsg,checkrun)
481                            return
482                        end
483                    end
484                else
485                    par_civ1.Mask=[];
486                end
487                mask=par_civ1.Mask;
488                maskoldname=maskname;
489            end
490        end
491        if strcmp(Param.ActionInput.ListCompareMode, 'PIV volume')
492            % Data.ListVarName=[Data.ListVarName 'Civ1_Z'];
493            % Data.Civ1_X=[];Data.Civ1_Y=[];Data.Civ1_Z=[];
494            % Data.Civ1_U=[];Data.Civ1_V=[];Data.Civ1_C=[];
495            % for ivol=1:NbSlice
496            %     % caluclate velocity data (y and v in indices, reverse to y component)
497            %     [xtable, ytable, utable, vtable, ctable, F, result_conv, errormsg] = civ (par_civ1);
498            %     if ~isempty(errormsg)
499            %         disp_uvmat('ERROR',errormsg,checkrun)
500            %         return
501            %     end
502            %     Data.Civ1_X=[Data.Civ1_X reshape(xtable,[],1)];
503            %     Data.Civ1_Y=[Data.Civ1_Y reshape(Param.Civ1.ImageHeight-ytable+1,[],1)];
504            %     Data.Civ1_Z=[Data.Civ1_Z ivol*ones(numel(xtable),1)];% z=image index in image coordinates
505            %     Data.Civ1_U=[Data.Civ1_U reshape(utable,[],1)];
506            %     Data.Civ1_V=[Data.Civ1_V reshape(-vtable,[],1)];
507            %     Data.Civ1_C=[Data.Civ1_C reshape(ctable,[],1)];
508            %     Data.Civ1_FF=[Data.Civ1_FF reshape(F,[],1)];
509            % end
510        else %usual PIV
511            % caluclate velocity data (y and v in indices, reverse to y component)
512            tstart_civ1=tic;
513            [xtable, ytable, utable, vtable, ctable, FF, result_conv, errormsg] = civ (par_civ1);
514            if ~isempty(errormsg)
515                disp_uvmat('ERROR',errormsg,checkrun)
516                return
517            end
518            Data.Civ1_X=reshape(xtable,[],1);
519            Data.Civ1_Y=reshape(par_civ1.ImageHeight-ytable+1,[],1);
520            Data.Civ1_U=reshape(utable,[],1);
521            Data.Civ1_V=reshape(-vtable,[],1);
522            Data.Civ1_C=reshape(ctable,[],1);
523            Data.Civ1_FF=reshape(FF,[],1);
524            time_civ1=toc(tstart_civ1);
525        end
526    else% we use existing Civ1 data
527        if exist('ncfile','var')
528            CivFile=ncfile;
529            [Data,tild,tild,errormsg]=nc2struct(CivFile,'ListGlobalAttribute','absolut_time_T0'); %look for the constant 'absolut_time_T0' to detect old civx data format
530            if ~isempty(errormsg)
531                disp_uvmat('ERROR',errormsg,checkrun)
532                return
533            end
534            [Data,tild,tild,errormsg]=nc2struct(CivFile);%read civ1 and fix1 data in the existing netcdf file
535        elseif isfield(Param,'Civ1_X')
536            Data.ListGlobalAttribute={};
537            Data.ListVarName={};
538            Data.VarDimName={};
539            Data.Civ1_X=Param.Civ1_X;
540            Data.Civ1_Y=Param.Civ1_Y;
541            Data.Civ1_U=Param.Civ1_U;
542            Data.Civ1_V=Param.Civ1_V;
543            Data.Civ1_FF=Param.Civ1_FF;
544        end
545    end
546   
547    %% Fix1
548    if isfield (Param.ActionInput,'Fix1')
549        disp('detect_false1 started')
550        if ~isfield (Param.ActionInput,'Civ1')% if we use existing Civ1, remove previous data beyond Civ1
551            Fix1_attr=find(strcmp('Fix1',Data.ListGlobalAttribute));
552            Data.ListGlobalAttribute(Fix1_attr)=[];
553            for ilist=1:numel(Fix1_attr)
554                Data=rmfield(Data,Data.ListGlobalAttribute{Fix1_attr(ilist)});
555            end
556        end
557        list_param=fieldnames(Param.ActionInput.Fix1)';
558        Fix1_param=regexprep(list_param,'^.+','Fix1_$0');% insert 'Fix1_' before  each string in ListFixParam
559        %indicate the values of all the global attributes in the output data
560        for ilist=1:length(list_param)
561            Data.(Fix1_param{ilist})=Param.ActionInput.Fix1.(list_param{ilist});
562        end
563        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Fix1_param];
564        Data.Civ1_FF=uint8(detect_false(Param.ActionInput.Fix1,Data.Civ1_C,Data.Civ1_U,Data.Civ1_V,Data.Civ1_FF));
565        Data.CivStage=2;
566    end
567    %% Patch1
568    if isfield (Param.ActionInput,'Patch1')
569        disp('patch1 started')
570         tstart_patch1=tic;
571       
572        % record the processing parameters of Patch1 as global attributes in the result nc file
573        list_param=fieldnames(Param.ActionInput.Patch1)';
574        list_param(strcmp('TestPatch1',list_param))=[];% remove 'TestPatch1' from the list of parameters
575        Patch1_param=regexprep(list_param,'^.+','Patch1_$0');% insert 'Patch1_' before  each parameter name
576        for ilist=1:length(list_param)
577            Data.(Patch1_param{ilist})=Param.ActionInput.Patch1.(list_param{ilist});
578        end
579        Data.CivStage=3;% record the new state of processing
580        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Patch1_param];
581       
582        % list the variables to record
583        nbvar=length(Data.ListVarName);
584        Data.ListVarName=[Data.ListVarName {'Civ1_U_smooth','Civ1_V_smooth','Civ1_SubRange','Civ1_NbCentres','Civ1_Coord_tps','Civ1_U_tps','Civ1_V_tps'}];
585        Data.VarDimName=[Data.VarDimName {'nb_vec_1','nb_vec_1',{'nb_coord','nb_bounds','nb_subdomain_1'},'nb_subdomain_1',...
586            {'nb_tps_1','nb_coord','nb_subdomain_1'},{'nb_tps_1','nb_subdomain_1'},{'nb_tps_1','nb_subdomain_1'}}];
587        Data.VarAttribute{nbvar+1}.Role='vector_x';
588        Data.VarAttribute{nbvar+2}.Role='vector_y';
589        Data.VarAttribute{nbvar+5}.Role='coord_tps';
590        Data.VarAttribute{nbvar+6}.Role='vector_x';
591        Data.VarAttribute{nbvar+7}.Role='vector_y';
592        Data.Civ1_U_smooth=Data.Civ1_U; % zeros(size(Data.Civ1_X));
593        Data.Civ1_V_smooth=Data.Civ1_V; %zeros(size(Data.Civ1_X));
594        if isfield(Data,'Civ1_FF')
595            ind_good=find(Data.Civ1_FF==0);
596        else
597            ind_good=1:numel(Data.Civ1_X);
598        end
599        if isempty(ind_good)
600                        disp_uvmat('ERROR','all vectors of civ1 are bad, check input parameters' ,checkrun)
601                        return
602        end
603
604        % perform Patch calculation using the UVMAT fct 'filter_tps'
605        [Data.Civ1_SubRange,Data.Civ1_NbCentres,Data.Civ1_Coord_tps,Data.Civ1_U_tps,Data.Civ1_V_tps,tild,Ures, Vres,tild,FFres]=...
606            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);
607        Data.Civ1_U_smooth(ind_good)=Ures;% take the interpolated (smoothed) velocity values for good vectors, keep civ1 data for the other
608        Data.Civ1_V_smooth(ind_good)=Vres;
609        Data.Civ1_FF(ind_good)=uint8(4*FFres);%set FF to value =4 for vectors eliminated by filter_tps
610        time_patch1=toc(tstart_patch1);
611        disp('patch1 performed')
612    end
613   
614    %% Civ2
615    if isfield (Param.ActionInput,'Civ2')
616        disp('civ2 started')
617        tstart_civ2=tic;
618        par_civ2=Param.ActionInput.Civ2;
619        if CheckInputFile % read input images (except in mode Test where it is introduced directly in Param.ActionInput.Civ1.ImageNameA and B)
620            par_civ2.ImageA=[];
621            par_civ2.ImageB=[];
622            if strcmp(Param.ActionInput.ListCompareMode,'displacement')
623                ImageName_A_Civ2=Param.ActionInput.RefFile;
624            else
625                ImageName_A_Civ2=fullfile_uvmat(RootPath_A,SubDir_A,RootFile_A,FileExt_A,NomType_A,i1_civ2,[],j1_civ2);
626            end
627            if strcmp(ImageName_A_Civ2,ImageName_A) && isequal(FrameIndex_A_Civ1(ifield),FrameIndex_A_Civ2(ifield))
628                par_civ2.ImageA=par_civ1.ImageA;
629            else
630                [par_civ2.ImageA,VideoObject_A] = read_image(ImageName_A_Civ2,FileType_A,VideoObject_A,FrameIndex_A_Civ2(ifield));
631            end
632            ImageName_B_Civ2=fullfile_uvmat(RootPath_B,SubDir_B,RootFile_B,FileExt_B,NomType_B,i2_civ2,[],j2_civ2);
633            if strcmp(ImageName_B_Civ2,ImageName_B) && isequal(FrameIndex_B_Civ1(ifield),FrameIndex_B_Civ2)
634                par_civ2.ImageB=par_civ1.ImageB;
635            else
636                [par_civ2.ImageB,VideoObject_B] = read_image(ImageName_B_Civ2,FileType_B,VideoObject_B,FrameIndex_B_Civ2(ifield));
637            end
638            [FileInfo_A,VideoObject_A]=get_file_info(ImageName_A_Civ2);
639            par_civ2.ImageWidth=FileInfo_A.Width;
640            par_civ2.ImageHeight=FileInfo_A.Height;
641            if isfield(par_civ2,'Grid')% grid points set as input file
642                if ischar(par_civ2.Grid)%read the grid file if the input is a file name
643                    par_civ2.Grid=dlmread(par_civ2.Grid);
644                    par_civ2.Grid(1,:)=[];%the first line must be removed (heading in the grid file)
645                end
646            else% automatic grid
647                minix=floor(par_civ2.Dx/2)-0.5;
648                maxix=minix+par_civ2.Dx*floor((par_civ2.ImageWidth-1)/par_civ2.Dx);
649                miniy=floor(par_civ2.Dy/2)-0.5;
650                maxiy=minix+par_civ2.Dy*floor((par_civ2.ImageHeight-1)/par_civ2.Dy);
651                [GridX,GridY]=meshgrid(minix:par_civ2.Dx:maxix,miniy:par_civ2.Dy:maxiy);
652                par_civ2.Grid(:,1)=reshape(GridX,[],1);
653                par_civ2.Grid(:,2)=reshape(GridY,[],1);
654            end
655        end
656       
657        % get the guess from patch1 or patch2 (case 'CheckCiv3')
658        if CheckInputFile % read input images (except in mode Test where it is introduced directly in Param.ActionInput.Civ1.ImageNameA and B)
659            if isfield (par_civ2,'CheckCiv3') && par_civ2.CheckCiv3 %get the guess from  patch2
660                SubRange= Data.Civ2_SubRange;
661                NbCentres=Data.Civ2_NbCentres;
662                Coord_tps=Data.Civ2_Coord_tps;
663                U_tps=Data.Civ2_U_tps;
664                V_tps=Data.Civ2_V_tps;
665                CivStage=Data.CivStage;%store the current CivStage
666                Civ1_Dt=Data.Civ2_Dt;
667                Data=[];%reinitialise the result structure Data
668                Data.ListGlobalAttribute={'Conventions','Program','CivStage'};
669                Data.Conventions='uvmat/civdata';% states the conventions used for the description of field variables and attributes
670                Data.Program='civ_series';
671                Data.CivStage=CivStage+1;%update the current civStage after reinitialisation of Data
672                Data.ListVarName={};
673                Data.VarDimName={};
674            else % get the guess from patch1
675                SubRange= Data.Civ1_SubRange;
676                NbCentres=Data.Civ1_NbCentres;
677                Coord_tps=Data.Civ1_Coord_tps;
678                U_tps=Data.Civ1_U_tps;
679                V_tps=Data.Civ1_V_tps;
680                Civ1_Dt=Data.Civ1_Dt;
681                Data.CivStage=4;
682            end
683        else
684            SubRange= par_civ2.Civ1_SubRange;
685            NbCentres=par_civ2.Civ1_NbCentres;
686            Coord_tps=par_civ2.Civ1_Coord_tps;
687            U_tps=par_civ2.Civ1_U_tps;
688            V_tps=par_civ2.Civ1_V_tps;
689            Civ1_Dt=par_civ2.Civ1_Dt;
690            Civ2_Dt=par_civ2.Civ1_Dt;
691            Data.ListVarName={};
692            Data.VarDimName={};
693        end
694        Shiftx=zeros(size(par_civ2.Grid,1),1);% shift expected from civ1 data
695        Shifty=zeros(size(par_civ2.Grid,1),1);
696        nbval=zeros(size(par_civ2.Grid,1),1);% nbre of interpolated values at each grid point (from the different patch subdomains)
697        if par_civ2.CheckDeformation
698            DUDX=zeros(size(par_civ2.Grid,1),1);
699            DUDY=zeros(size(par_civ2.Grid,1),1);
700            DVDX=zeros(size(par_civ2.Grid,1),1);
701            DVDY=zeros(size(par_civ2.Grid,1),1);
702        end
703        NbSubDomain=size(SubRange,3);
704        for isub=1:NbSubDomain% for each sub-domain of Patch1
705            nbvec_sub=NbCentres(isub);% nbre of Civ vectors in the subdomain
706            ind_sel=find(par_civ2.Grid(:,1)>=SubRange(1,1,isub) & par_civ2.Grid(:,1)<=SubRange(1,2,isub) &...
707                par_civ2.Grid(:,2)>=SubRange(2,1,isub) & par_civ2.Grid(:,2)<=SubRange(2,2,isub));% grid points in the subdomain
708            if ~isempty(ind_sel)
709                epoints = par_civ2.Grid(ind_sel,:);% coordinates of interpolation sites (measurement grids)
710                ctrs=Coord_tps(1:nbvec_sub,:,isub) ;%(=initial points) ctrs
711                nbval(ind_sel)=nbval(ind_sel)+1;% records the number of values for each interpolation point (in case of subdomain overlap)
712                EM = tps_eval(epoints,ctrs);% thin plate spline (tps) coefficient
713                Shiftx(ind_sel)=Shiftx(ind_sel)+EM*U_tps(1:nbvec_sub+3,isub);%velocity shift estimated by tps from civ1
714                Shifty(ind_sel)=Shifty(ind_sel)+EM*V_tps(1:nbvec_sub+3,isub);
715                if par_civ2.CheckDeformation
716                    [EMDX,EMDY] = tps_eval_dxy(epoints,ctrs);%2D matrix of distances between extrapolation points epoints and spline centres (=site points) ctrs
717                    DUDX(ind_sel)=DUDX(ind_sel)+EMDX*U_tps(1:nbvec_sub+3,isub);
718                    DUDY(ind_sel)=DUDY(ind_sel)+EMDY*U_tps(1:nbvec_sub+3,isub);
719                    DVDX(ind_sel)=DVDX(ind_sel)+EMDX*V_tps(1:nbvec_sub+3,isub);
720                    DVDY(ind_sel)=DVDY(ind_sel)+EMDY*V_tps(1:nbvec_sub+3,isub);
721                end
722            end
723        end
724        i1=i1_series_Civ2(ifield);
725        i2=i1;
726        if ~isempty(i2_series_Civ2)
727            i2=i2_series_Civ1(ifield);
728        end
729        j1=1;
730        if ~isempty(j1_series_Civ2)
731            j1=j1_series_Civ1(ifield);
732        end
733        j2=j1;
734        if ~isempty(j2_series_Civ2)
735            j2=j2_series_Civ2(ifield);
736        end
737        if par_civ2.CheckMask&&~isempty(par_civ2.Mask)
738            if isfield(par_civ2,'NbSlice')
739                [RootPath_mask,SubDir_mask,RootFile_mask,~,~,~,~,Ext_mask]=fileparts_uvmat(Param.ActionInput.Civ2.Mask);
740                i1_mask=mod(i1-1,par_civ2.NbSlice)+1;
741                maskname=fullfile_uvmat(RootPath_mask,SubDir_mask,RootFile_mask,Ext_mask,'_1',i1_mask);
742            else
743                maskname=Param.ActionInput.Civ2.Mask;
744            end
745            if strcmp(maskoldname,maskname)% mask exist, not already read in civ1
746                par_civ2.Mask=mask; %use mask already opened
747            else
748                if exist(maskname,'file')
749                try
750                    par_civ2.Mask=imread(maskname);%update the mask, an store it for future use
751                catch ME
752                    if ~isempty(ME.message)
753                        errormsg=['error reading input image: ' ME.message];
754                        disp_uvmat('ERROR',errormsg,checkrun)
755                        return
756                    end
757                end
758                else
759                    par_civ2.Mask=[];
760                end
761                mask=par_civ2.Mask;
762                maskoldname=maskname;
763            end
764        end
765
766        if CheckInputFile % else Dt given by par_civ2
767            if strcmp(Param.ActionInput.ListCompareMode,'displacement')
768                Civ1_Dt=1;
769                Civ2_Dt=1;
770            else
771                Civ2_Dt=Time(i2_civ2+1,j2_civ2+1)-Time(i1_civ2+1,j1_civ2+1);
772            end
773        end
774        par_civ2.SearchBoxShift=(Civ2_Dt/Civ1_Dt)*[Shiftx(nbval>=1)./nbval(nbval>=1) Shifty(nbval>=1)./nbval(nbval>=1)];
775        % shift the grid points by half the expected shift to provide the correlation box position in image A
776        par_civ2.Grid=[par_civ2.Grid(nbval>=1,1)-par_civ2.SearchBoxShift(:,1)/2 par_civ2.Grid(nbval>=1,2)-par_civ2.SearchBoxShift(:,2)/2];
777        if par_civ2.CheckDeformation
778            par_civ2.DUDX=DUDX(nbval>=1)./nbval(nbval>=1);
779            par_civ2.DUDY=DUDY(nbval>=1)./nbval(nbval>=1);
780            par_civ2.DVDX=DVDX(nbval>=1)./nbval(nbval>=1);
781            par_civ2.DVDY=DVDY(nbval>=1)./nbval(nbval>=1);
782        end
783       
784        % calculate velocity data (y and v in image indices, reverse to y component)
785       
786        [xtable, ytable, utable, vtable, ctable, FF,result_conv,errormsg] = civ (par_civ2);
787       
788        list_param=(fieldnames(Param.ActionInput.Civ2))';
789        list_param(strcmp('TestCiv2',list_param))=[];% remove the parameter TestCiv2 from the list
790        Civ2_param=regexprep(list_param,'^.+','Civ2_$0');% insert 'Civ2_' before  each string in list_param
791        Civ2_param=[{'Civ2_ImageA','Civ2_ImageB','Civ2_Time','Civ2_Dt'} Civ2_param]; %insert the names of the two input images
792        %indicate the values of all the global attributes in the output data
793        if exist('ImageName_A','var')
794            Data.Civ2_ImageA=ImageName_A;
795            Data.Civ2_ImageB=ImageName_B;
796            if strcmp(Param.ActionInput.ListCompareMode,'displacement')
797                Data.Civ2_Time=Time(i2_civ2+1,j2_civ2+1);% the Time is the Time of the secodn image
798                Data.Civ2_Dt=1;% Time interval is 1, to yield displacement instead of velocity=displacement/Dt at reading
799            else
800                Data.Civ2_Time=(Time(i2_civ2+1,j2_civ2+1)+Time(i1_civ2+1,j1_civ2+1))/2;
801                Data.Civ2_Dt=Civ2_Dt;
802            end
803        end
804        for ilist=1:length(list_param)
805            Data.(Civ2_param{4+ilist})=Param.ActionInput.Civ2.(list_param{ilist});
806        end
807        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Civ2_param];
808       
809        nbvar=numel(Data.ListVarName);
810        % define the Civ2 variable (if Civ2 data are not replaced from previous calculation)
811        if isempty(find(strcmp('Civ2_X',Data.ListVarName),1))
812            Data.ListVarName=[Data.ListVarName {'Civ2_X','Civ2_Y','Civ2_U','Civ2_V','Civ2_C','Civ2_FF'}];%  cell array containing the names of the fields to record
813            Data.VarDimName=[Data.VarDimName {'nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2','nb_vec_2'}];
814            Data.VarAttribute{nbvar+1}.Role='coord_x';
815            Data.VarAttribute{nbvar+2}.Role='coord_y';
816            Data.VarAttribute{nbvar+3}.Role='vector_x';
817            Data.VarAttribute{nbvar+4}.Role='vector_y';
818            Data.VarAttribute{nbvar+5}.Role='ancillary';
819            Data.VarAttribute{nbvar+6}.Role='errorflag';
820        end
821        Data.Civ2_X=reshape(xtable,[],1);
822        Data.Civ2_Y=reshape(size(par_civ2.ImageA,1)-ytable+1,[],1);
823        Data.Civ2_U=reshape(utable,[],1);
824        Data.Civ2_V=reshape(-vtable,[],1);
825        Data.Civ2_C=reshape(ctable,[],1);
826        Data.Civ2_FF=reshape(FF,[],1);
827        disp('civ2 performed')
828        time_civ2=toc(tstart_civ2);
829    elseif ~isfield(Data,'ListVarName') % we start there, using existing Civ2 data
830        if exist('ncfile','var')
831            CivFile=ncfile;
832            [Data,tild,tild,errormsg]=nc2struct(CivFile);%read civ1 and detect_false1 data in the existing netcdf file
833            if ~isempty(errormsg)
834                disp_uvmat('ERROR',errormsg,checkrun)
835                return
836            end
837        end
838    end
839   
840    %% Fix2
841    if isfield (Param.ActionInput,'Fix2')
842        disp('detect_false2 started')
843        list_param=fieldnames(Param.ActionInput.Fix2)';
844        Fix2_param=regexprep(list_param,'^.+','Fix2_$0');% insert 'Fix1_' before  each string in ListFixParam
845        %indicate the values of all the global attributes in the output data
846        for ilist=1:length(list_param)
847            Data.(Fix2_param{ilist})=Param.ActionInput.Fix2.(list_param{ilist});
848        end
849        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Fix2_param];     
850        Data.Civ2_FF=double(detect_false(Param.ActionInput.Fix2,Data.Civ2_C,Data.Civ2_U,Data.Civ2_V,Data.Civ2_FF));
851        Data.CivStage=Data.CivStage+1;
852    end
853   
854    %% Patch2
855    if isfield (Param.ActionInput,'Patch2')
856        disp('patch2 started')
857        tstart_patch2=tic;
858        list_param=fieldnames(Param.ActionInput.Patch2)';
859        list_param(strcmp('TestPatch2',list_param))=[];% remove the parameter TestCiv1 from the list
860        Patch2_param=regexprep(list_param,'^.+','Patch2_$0');% insert 'Fix1_' before  each string in ListFixParam
861        %indicate the values of all the global attributes in the output data
862        for ilist=1:length(list_param)
863            Data.(Patch2_param{ilist})=Param.ActionInput.Patch2.(list_param{ilist});
864        end
865        Data.ListGlobalAttribute=[Data.ListGlobalAttribute Patch2_param];
866       
867        nbvar=length(Data.ListVarName);
868        Data.ListVarName=[Data.ListVarName {'Civ2_U_smooth','Civ2_V_smooth','Civ2_SubRange','Civ2_NbCentres','Civ2_Coord_tps','Civ2_U_tps','Civ2_V_tps'}];
869        Data.VarDimName=[Data.VarDimName {'nb_vec_2','nb_vec_2',{'nb_coord','nb_bounds','nb_subdomain_2'},{'nb_subdomain_2'},...
870            {'nb_tps_2','nb_coord','nb_subdomain_2'},{'nb_tps_2','nb_subdomain_2'},{'nb_tps_2','nb_subdomain_2'}}];
871       
872        Data.VarAttribute{nbvar+1}.Role='vector_x';
873        Data.VarAttribute{nbvar+2}.Role='vector_y';
874        Data.VarAttribute{nbvar+5}.Role='coord_tps';
875        Data.VarAttribute{nbvar+6}.Role='vector_x';
876        Data.VarAttribute{nbvar+7}.Role='vector_y';
877        Data.Civ2_U_smooth=Data.Civ2_U;
878        Data.Civ2_V_smooth=Data.Civ2_V;
879        if isfield(Data,'Civ2_FF')
880            ind_good=find(Data.Civ2_FF==0);
881        else
882            ind_good=1:numel(Data.Civ2_X);
883        end
884                if isempty(ind_good)
885                        disp_uvmat('ERROR','all vectors of civ2 are bad, check input parameters' ,checkrun)
886                        return
887                end
888             
889        [Data.Civ2_SubRange,Data.Civ2_NbCentres,Data.Civ2_Coord_tps,Data.Civ2_U_tps,Data.Civ2_V_tps,tild,Ures,Vres,tild,FFres]=...
890            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);
891        Data.Civ2_U_smooth(ind_good)=Ures;
892        Data.Civ2_V_smooth(ind_good)=Vres;
893        Data.Civ2_FF(ind_good)=uint8(4*FFres);
894        Data.CivStage=Data.CivStage+1;
895        time_patch2=toc(tstart_patch2);
896        disp('patch2 performed')
897    end
898   
899    %% write result in a netcdf file if requested
900    if CheckOutputFile
901        errormsg=struct2nc(ncfile_out,Data);
902        if isempty(errormsg)
903            disp([ncfile_out ' written'])
904            %[success,msg] = fileattrib(ncfile_out ,'+w','g');% done in struct2nc
905        else
906            disp(errormsg)
907        end
908        time_total=toc(tstart);
909        disp(['ellapsed time ' num2str(time_total/60,2) ' minutes'])
910        disp(['time civ1 ' num2str(time_civ1,2) ' s'])
911        disp(['time patch1 ' num2str(time_patch1,2) ' s'])
912        disp(['time civ2 ' num2str(time_civ2,2) ' s'])
913        disp(['time patch2 ' num2str(time_patch2,2) ' s'])
914        if exist('time_input','var')
915            disp(['time image reading ' num2str(time_input,2) ' s'])
916            disp(['time other ' num2str((time_total-time_input-time_civ1-time_patch1-time_civ2-time_patch2),2) ' s'])
917        end
918    end
919end
920
921
922% 'civ': function piv.m adapted from PIVlab http://pivlab.blogspot.com/
923%--------------------------------------------------------------------------
924% function [xtable ytable utable vtable typevector] = civ (image1,image2,ibx,iby step, subpixfinder, mask, roi)
925%
926% OUTPUT:
927% xtable: set of x coordinates
928% ytable: set of y coordinates
929% utable: set of u displacements (along x)
930% vtable: set of v displacements (along y)
931% ctable: max image correlation for each vector
932% typevector: set of flags, =1 for good, =0 for NaN vectors
933%
934%INPUT:
935% par_civ: structure of input parameters, with fields:
936%  .ImageA: first image for correlation (matrix)
937%  .ImageB: second image for correlation(matrix)
938%  .CorrBoxSize: 1,2 vector giving the size of the correlation box in x and y
939%  .SearchBoxSize:  1,2 vector giving the size of the search box in x and y
940%  .SearchBoxShift: 1,2 vector or 2 column matrix (for civ2) giving the shift of the search box in x and y
941%  .CorrSmooth: =1 or 2 determines the choice of the sub-pixel determination of the correlation max
942%  .ImageWidth: nb of pixels of the image in x
943%  .Dx, Dy: mesh for the PIV calculation
944%  .Grid: grid giving the PIV calculation points (alternative to .Dx .Dy): centres of the correlation boxes in Image A
945%  .Mask: name of a mask file or mask image matrix itself
946%  .MinIma: thresholds for image luminosity
947%  .MaxIma
948%  .CheckDeformation=1 for subpixel interpolation and image deformation (linear transform)
949%  .DUDX: matrix of deformation obtained from patch at each grid point
950%  .DUDY
951%  .DVDX:
952%  .DVDY
953
954function [xtable,ytable,utable,vtable,ctable,F,result_conv,errormsg] = civ (par_civ)
955
956%% prepare measurement grid
957if isfield(par_civ,'Grid')% grid points set as input, central positions of the sub-images in image A
958    if ischar(par_civ.Grid)%read the grid file if the input is a file name (grid in x, y image coordinates)
959        par_civ.Grid=dlmread(par_civ.Grid);
960        par_civ.Grid(1,:)=[];%the first line must be removed (heading in the grid file)
961    end
962    % else par_civ.Grid is already an array, no action here
963else% automatic grid in x, y image coordinates
964    minix=floor(par_civ.Dx/2)-0.5;
965    maxix=minix+par_civ.Dx*floor((par_civ.ImageWidth-1)/par_civ.Dx);
966    miniy=floor(par_civ.Dy/2)-0.5;% first automatic grid point at half the mesh Dy
967    maxiy=minix+par_civ.Dy*floor((par_civ.ImageHeight-1)/par_civ.Dy);
968    [GridX,GridY]=meshgrid(minix:par_civ.Dx:maxix,miniy:par_civ.Dy:maxiy);
969    par_civ.Grid(:,1)=reshape(GridX,[],1);
970    par_civ.Grid(:,2)=reshape(GridY,[],1);% increases with array index
971end
972nbvec=size(par_civ.Grid,1);
973
974%% prepare correlation and search boxes
975ibx2=floor(par_civ.CorrBoxSize(1)/2);
976iby2=floor(par_civ.CorrBoxSize(2)/2);
977isx2=floor(par_civ.SearchBoxSize(1)/2);
978isy2=floor(par_civ.SearchBoxSize(2)/2);
979shiftx=round(par_civ.SearchBoxShift(:,1));%use the input shift estimate, rounded to the next integer value
980shifty=-round(par_civ.SearchBoxShift(:,2));% sign minus because image j index increases when y decreases
981if numel(shiftx)==1% case of a unique shift for the whole field( civ1)
982    shiftx=shiftx*ones(nbvec,1);
983    shifty=shifty*ones(nbvec,1);
984end
985
986%% Array initialisation and default output  if par_civ.CorrSmooth=0 (just the grid calculated, no civ computation)
987xtable=round(par_civ.Grid(:,1)+0.5)-0.5;
988ytable=round(par_civ.ImageHeight-par_civ.Grid(:,2)+0.5)-0.5;% y index corresponding to the position in image coordiantes
989utable=shiftx;%zeros(nbvec,1);
990vtable=shifty;%zeros(nbvec,1);
991ctable=zeros(nbvec,1);
992F=zeros(nbvec,1);
993result_conv=[];
994errormsg='';
995
996%% prepare mask
997if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
998    if strcmp(par_civ.Mask,'all')
999        return    % get the grid only, no civ calculation
1000    elseif ischar(par_civ.Mask)
1001        par_civ.Mask=imread(par_civ.Mask);% read the mask if not allready done
1002    end
1003end
1004check_MinIma=isfield(par_civ,'MinIma');% test for image luminosity threshold
1005check_MaxIma=isfield(par_civ,'MaxIma') && ~isempty(par_civ.MaxIma);
1006
1007par_civ.ImageA=sum(double(par_civ.ImageA),3);%sum over rgb component for color images
1008par_civ.ImageB=sum(double(par_civ.ImageB),3);
1009[npy_ima npx_ima]=size(par_civ.ImageA);
1010if ~isequal(size(par_civ.ImageB),[npy_ima npx_ima])
1011    errormsg='image pair with unequal size';
1012    return
1013end
1014
1015%% Apply mask
1016% Convention for mask, IDEAS NOT IMPLEMENTED
1017% mask >200 : velocity calculated
1018%  200 >=mask>150;velocity not calculated, interpolation allowed (bad spots)
1019% 150>=mask >100: velocity not calculated, nor interpolated
1020%  100>=mask> 20: velocity not calculated, impermeable (no flux through mask boundaries)
1021%  20>=mask: velocity=0
1022checkmask=0;
1023MinA=min(min(par_civ.ImageA));
1024%MinB=min(min(par_civ.ImageB));
1025%check_undefined=false(size(par_civ.ImageA));
1026if isfield(par_civ,'Mask') && ~isempty(par_civ.Mask)
1027    checkmask=1;
1028    if ~isequal(size(par_civ.Mask),[npy_ima npx_ima])
1029        errormsg='mask must be an image with the same size as the images';
1030        return
1031    end
1032    check_undefined=(par_civ.Mask<200 & par_civ.Mask>=20 );
1033end
1034
1035%% compute image correlations: MAINLOOP on velocity vectors
1036corrmax=0;
1037sum_square=1;% default
1038mesh=1;% default
1039CheckDeformation=isfield(par_civ,'CheckDeformation')&& par_civ.CheckDeformation==1;
1040if CheckDeformation
1041    mesh=0.25;%mesh in pixels for subpixel image interpolation (x 4 in each direction)
1042    par_civ.CorrSmooth=2;% use SUBPIX2DGAUSS (take into account more points near the max)
1043end
1044
1045if par_civ.CorrSmooth~=0 % par_civ.CorrSmooth=0 implies no civ computation (just input image and grid points given)
1046    for ivec=1:nbvec
1047        iref=round(par_civ.Grid(ivec,1)+0.5);% xindex on the image A for the middle of the correlation box
1048        jref=round(par_civ.ImageHeight-par_civ.Grid(ivec,2)+0.5);%  j index  for the middle of the correlation box in the image A
1049        F(ivec)=0;
1050        subrange1_x=iref-ibx2:iref+ibx2;% x indices defining the first subimage
1051        subrange1_y=jref-iby2:jref+iby2;% y indices defining the first subimage
1052        subrange2_x=iref+shiftx(ivec)-isx2:iref+shiftx(ivec)+isx2;%x indices defining the second subimage
1053        subrange2_y=jref+shifty(ivec)-isy2:jref+shifty(ivec)+isy2;%y indices defining the second subimage
1054        image1_crop=MinA*ones(numel(subrange1_y),numel(subrange1_x));% default value=min of image A
1055        image2_crop=MinA*ones(numel(subrange2_y),numel(subrange2_x));% default value=min of image A
1056        check1_x=subrange1_x>=1 & subrange1_x<=par_civ.ImageWidth;% check which points in the subimage 1 are contained in the initial image 1
1057        check1_y=subrange1_y>=1 & subrange1_y<=par_civ.ImageHeight;
1058        check2_x=subrange2_x>=1 & subrange2_x<=par_civ.ImageWidth;% check which points in the subimage 2 are contained in the initial image 2
1059        check2_y=subrange2_y>=1 & subrange2_y<=par_civ.ImageHeight;
1060        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
1061        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
1062        if checkmask
1063            mask1_crop=ones(numel(subrange1_y),numel(subrange1_x));% default value=1 for mask
1064            mask2_crop=ones(numel(subrange2_y),numel(subrange2_x));% default value=1 for mask
1065            mask1_crop(check1_y,check1_x)=check_undefined(subrange1_y(check1_y),subrange1_x(check1_x));%extract a mask subimage (correlation box) from image A
1066            mask2_crop(check2_y,check2_x)=check_undefined(subrange2_y(check2_y),subrange2_x(check2_x));%extract a mask subimage (search box) from image B
1067            sizemask=sum(sum(mask1_crop))/(numel(subrange1_y)*numel(subrange1_x));%size of the masked part relative to the correlation sub-image
1068            if sizemask > 1/2% eliminate point if more than half of the correlation box is masked
1069                F(ivec)=1; %
1070                utable(ivec)=NaN;
1071                vtable(ivec)=NaN;
1072            else
1073                image1_crop=image1_crop.*~mask1_crop;% put to zero the masked pixels (mask1_crop='true'=1)
1074                image2_crop=image2_crop.*~mask2_crop;
1075                image1_mean=mean(mean(image1_crop))/(1-sizemask);
1076                image2_mean=mean(mean(image2_crop))/(1-sizemask);
1077            end
1078        else
1079            image1_mean=mean(mean(image1_crop));
1080            image2_mean=mean(mean(image2_crop));
1081        end
1082        %threshold on image minimum
1083        if F(ivec)~=1
1084            if check_MinIma && (image1_mean < par_civ.MinIma || image2_mean < par_civ.MinIma)
1085                F(ivec)=1;
1086                %threshold on image maximum
1087            elseif check_MaxIma && (image1_mean > par_civ.MaxIma || image2_mean > par_civ.MaxIma)
1088                F(ivec)=1;
1089            end
1090            if F(ivec)==1
1091                utable(ivec)=NaN;
1092                vtable(ivec)=NaN;
1093            else
1094                %mask
1095                if checkmask
1096                    image1_crop=(image1_crop-image1_mean).*~mask1_crop;%substract the mean, put to zero the masked parts
1097                    image2_crop=(image2_crop-image2_mean).*~mask2_crop;
1098                else
1099                    image1_crop=(image1_crop-image1_mean);
1100                    image2_crop=(image2_crop-image2_mean);
1101                end
1102                %deformation
1103                if CheckDeformation
1104                    xi=(1:mesh:size(image1_crop,2));
1105                    yi=(1:mesh:size(image1_crop,1))';
1106                    [XI,YI]=meshgrid(xi-ceil(size(image1_crop,2)/2),yi-ceil(size(image1_crop,1)/2));
1107                    XIant=XI-par_civ.DUDX(ivec)*XI+par_civ.DUDY(ivec)*YI+ceil(size(image1_crop,2)/2);
1108                    YIant=YI+par_civ.DVDX(ivec)*XI-par_civ.DVDY(ivec)*YI+ceil(size(image1_crop,1)/2);
1109                    image1_crop=interp2(image1_crop,XIant,YIant);
1110                    image1_crop(isnan(image1_crop))=0;
1111                    xi=(1:mesh:size(image2_crop,2));
1112                    yi=(1:mesh:size(image2_crop,1))';
1113                    image2_crop=interp2(image2_crop,xi,yi,'*spline');
1114                    image2_crop(isnan(image2_crop))=0;
1115                end
1116                sum_square=sum(sum(image1_crop.*image1_crop));
1117                %reference: Oliver Pust, PIV: Direct Cross-Correlation
1118                result_conv= conv2(image2_crop,flipdim(flipdim(image1_crop,2),1),'valid');
1119                corrmax= max(max(result_conv));
1120                result_conv=(result_conv/corrmax)*255; %normalize, peak=always 255
1121                %Find the correlation max, at 255
1122                [y,x] = find(result_conv==255,1);
1123                subimage2_crop=image2_crop(y:y+2*iby2/mesh,x:x+2*ibx2/mesh);%subimage of image 2 corresponding to the optimum displacement of first image
1124                sum_square=sum_square*sum(sum(subimage2_crop.*subimage2_crop));% product of variances of image 1 and 2
1125                sum_square=sqrt(sum_square);% srt of the variance product to normalise correlation
1126                if ~isempty(y) && ~isempty(x)
1127                    try
1128                        if par_civ.CorrSmooth==1
1129                            [vector,F(ivec)] = SUBPIXGAUSS (result_conv,x,y);
1130                        elseif par_civ.CorrSmooth==2
1131                            [vector,F(ivec)] = SUBPIX2DGAUSS (result_conv,x,y);
1132                        else
1133                            [vector,F(ivec)] = quadr_fit(result_conv,x,y);
1134                        end
1135                        utable(ivec)=vector(1)*mesh+shiftx(ivec);
1136                        vtable(ivec)=vector(2)*mesh+shifty(ivec);
1137                        xtable(ivec)=iref+utable(ivec)/2-0.5;% convec flow (velocity taken at the point middle from imgae 1 and 2)
1138                        ytable(ivec)=jref+vtable(ivec)/2-0.5;% and position of pixel 1=0.5 (convention for image coordinates=0 at the edge)
1139                        iref=round(xtable(ivec)+0.5);% nearest image index for the middle of the vector
1140                        jref=round(ytable(ivec)+0.5);
1141                        % eliminate vectors located in the mask
1142                        if  checkmask && (iref<1 || jref<1 ||iref>npx_ima || jref>npy_ima ||( par_civ.Mask(jref,iref)<200 && par_civ.Mask(jref,iref)>=100))
1143                            utable(ivec)=0;
1144                            vtable(ivec)=0;
1145                            F(ivec)=1;
1146                        end
1147                        ctable(ivec)=corrmax/sum_square;% correlation value
1148                    catch ME
1149                        F(ivec)=1;
1150                        disp(ME.message)
1151                    end
1152                else
1153                    F(ivec)=1;
1154                end
1155            end
1156        end
1157    end
1158end
1159result_conv=result_conv*corrmax/(255*sum_square);% keep the last correlation matrix for output
1160
1161%------------------------------------------------------------------------
1162% --- Find the maximum of the correlation function after interpolation
1163% OUPUT:
1164% vector = optimum displacement vector with subpixel correction
1165% F =flag: =0 OK
1166%           =-2 , warning: max too close to the edge of the search box (1 pixel margin)
1167% INPUT:
1168% x,y: position of the maximum correlation at integer values
1169
1170function [vector,F] = SUBPIXGAUSS (result_conv,x,y)
1171%------------------------------------------------------------------------
1172% vector=[0 0]; %default
1173F=0;
1174[npy,npx]=size(result_conv);
1175result_conv(result_conv<1)=1; %set to 1 correlation values smaller than 1  (=0 by discretisation, to avoid divergence in the log)
1176%the following 8 lines are copyright (c) 1998, Uri Shavit, Roi Gurka, Alex Liberzon, Technion ??? Israel Institute of Technology
1177%http://urapiv.wordpress.com
1178peaky = y;
1179if y < npy && y > 1
1180    f0 = log(result_conv(y,x));
1181    f1 = log(result_conv(y-1,x));
1182    f2 = log(result_conv(y+1,x));
1183    peaky = peaky+ (f1-f2)/(2*f1-4*f0+2*f2);
1184else
1185    F=1; % warning flag for vector truncated by the limited search box
1186end
1187peakx=x;
1188if x < npx-1 && x > 1
1189    f0 = log(result_conv(y,x));
1190    f1 = log(result_conv(y,x-1));
1191    f2 = log(result_conv(y,x+1));
1192    peakx = peakx+ (f1-f2)/(2*f1-4*f0+2*f2);
1193else
1194    F=1; % warning flag for vector truncated by the limited search box
1195end
1196vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
1197
1198%------------------------------------------------------------------------
1199% --- Find the maximum of the correlation function after interpolation
1200function [vector,F] = SUBPIX2DGAUSS (result_conv,x,y)
1201%------------------------------------------------------------------------
1202% vector=[0 0]; %default
1203F=1;
1204peaky=y;
1205peakx=x;
1206result_conv(result_conv<1)=1; %set to 1 correlation values smaller than 1 (to avoid divergence in the log)
1207[npy,npx]=size(result_conv);
1208if (x < npx) && (y < npy) && (x > 1) && (y > 1)
1209    F=0;
1210    for i=-1:1
1211        for j=-1:1
1212            %following 15 lines based on
1213            %H. Nobach ??? M. Honkanen (2005)
1214            %Two-dimensional Gaussian regression for sub-pixel displacement
1215            %estimation in particle image velocimetry or particle position
1216            %estimation in particle tracking velocimetry
1217            %Experiments in Fluids (2005) 38: 511???515
1218            c10(j+2,i+2)=i*log(result_conv(y+j, x+i));
1219            c01(j+2,i+2)=j*log(result_conv(y+j, x+i));
1220            c11(j+2,i+2)=i*j*log(result_conv(y+j, x+i));
1221            c20(j+2,i+2)=(3*i^2-2)*log(result_conv(y+j, x+i));
1222            c02(j+2,i+2)=(3*j^2-2)*log(result_conv(y+j, x+i));
1223        end
1224    end
1225    c10=(1/6)*sum(sum(c10));
1226    c01=(1/6)*sum(sum(c01));
1227    c11=(1/4)*sum(sum(c11));
1228    c20=(1/6)*sum(sum(c20));
1229    c02=(1/6)*sum(sum(c02));
1230    deltax=(c11*c01-2*c10*c02)/(4*c20*c02-c11^2);
1231    deltay=(c11*c10-2*c01*c20)/(4*c20*c02-c11^2);
1232    if abs(deltax)<1
1233        peakx=x+deltax;
1234    end
1235    if abs(deltay)<1
1236        peaky=y+deltay;
1237    end
1238end
1239vector=[peakx-floor(npx/2)-1 peaky-floor(npy/2)-1];
1240
1241%------------------------------------------------------------------------
1242% --- Find the maximum of the correlation function after quadratic interpolation
1243function [vector,F] = quadr_fit(result_conv,x,y)
1244[npy,npx]=size(result_conv);
1245if x<4 || y<4 || npx-x<4 ||npy-y <4
1246    F=1;
1247    vector=[x y];
1248else
1249    F=0;
1250    x_ind=x-4:x+4;
1251    y_ind=y-4:y+4;
1252    x_vec=0.25*(x_ind-x);
1253    y_vec=0.25*(y_ind-y);
1254    [X,Y]=meshgrid(x_vec,y_vec);
1255    coord=[reshape(X,[],1) reshape(Y,[],1)];
1256    result_conv=reshape(result_conv(y_ind,x_ind),[],1);
1257   
1258   
1259    % n=numel(X);
1260    % x=[X Y];
1261    % X=X-0.5;
1262    % Y=Y+0.5;
1263    % y = (X.*X+2*Y.*Y+X.*Y+6) + 0.1*rand(n,1);
1264    p = polyfitn(coord,result_conv,2);
1265    A(1,1)=2*p.Coefficients(1);
1266    A(1,2)=p.Coefficients(2);
1267    A(2,1)=p.Coefficients(2);
1268    A(2,2)=2*p.Coefficients(4);
1269    vector=[x y]'-A\[p.Coefficients(3) p.Coefficients(5)]';
1270    vector=vector'-[floor(npx/2) floor(npy/2)]-1 ;
1271    % zg = polyvaln(p,coord);
1272    % figure
1273    % surf(x_vec,y_vec,reshape(zg,9,9))
1274    % hold on
1275    % plot3(X,Y,reshape(result_conv,9,9),'o')
1276    % hold off
1277end
1278
1279
1280function FF=detect_false(Param,C,U,V,FFIn)
1281FF=FFIn;%default, good vectors
1282% FF=1, for correlation max at edge, not set in this function
1283% FF=2, for too small correlation
1284% FF=3, for velocity outside bounds
1285% FF=4 for exclusion by difference with the smoothed field, set by call to function filter_tps
1286
1287if isfield (Param,'MinCorr')
1288     FF(C<Param.MinCorr & FFIn==0)=2;
1289end
1290if (isfield(Param,'MinVel')&&~isempty(Param.MinVel))||(isfield (Param,'MaxVel')&&~isempty(Param.MaxVel))
1291    Umod= U.*U+V.*V;
1292    if isfield (Param,'MinVel')&&~isempty(Param.MinVel)
1293        U2Min=Param.MinVel*Param.MinVel;
1294        FF(Umod<U2Min & FFIn==0)=3;
1295    end
1296    if isfield (Param,'MaxVel')&&~isempty(Param.MaxVel)
1297         U2Max=Param.MaxVel*Param.MaxVel;
1298        FF(Umod>U2Max & FFIn==0)=3;
1299    end
1300end
1301
1302%------------------------------------------------------------------------
1303% --- determine the list of index pairs of processing file
1304function [i1_series,i2_series,j1_series,j2_series,check_bounds,NomTypeNc]=...
1305    find_pair_indices(str_civ,i_series,j_series,MinIndex_i,MaxIndex_i,MinIndex_j,MaxIndex_j)
1306%------------------------------------------------------------------------
1307i1_series=i_series;% set of first image indexes
1308i2_series=i_series;
1309j1_series=j_series;%ones(size(i_series));% set of first image numbers
1310j2_series=j_series;%ones(size(i_series));
1311r=regexp(str_civ,'^\D(?<ind>[i|j])=( -| )(?<num1>\d+)\|(?<num2>\d+)','names');
1312if ~isempty(r)
1313    mode=['D' r.ind];
1314    ind1=str2num(r.num1);
1315    ind2=str2num(r.num2);
1316else
1317    mode='j1-j2';
1318    r=regexp(str_civ,'^j= (?<num1>[a-z])-(?<num2>[a-z])','names');
1319    if ~isempty(r)
1320        NomTypeNc='_1ab';
1321    else
1322        r=regexp(str_civ,'^j= (?<num1>[A-Z])-(?<num2>[A-Z])','names');
1323        if ~isempty(r)
1324            NomTypeNc='_1AB';
1325        else
1326            r=regexp(str_civ,'^j= (?<num1>\d+)-(?<num2>\d+)','names');
1327            if ~isempty(r)
1328                NomTypeNc='_1_1-2';
1329            end
1330        end
1331    end
1332    if isempty(r)
1333        display('wrong pair mode input option')
1334    else
1335        ind1=stra2num(r.num1);
1336        ind2=stra2num(r.num2);
1337    end
1338end
1339switch mode
1340    case 'Di'
1341        i1_series=i_series-ind1;% set of first image numbers
1342        i2_series=i_series+ind2;
1343        check_bounds=i1_series<MinIndex_i | i2_series>MaxIndex_i;
1344        if isempty(j_series)
1345            NomTypeNc='_1-2';
1346        else
1347            j1_series=j_series;
1348            j2_series=j_series;
1349            NomTypeNc='_1-2_1';
1350        end
1351    case 'Dj'
1352        j1_series=j_series-ind1;
1353        j2_series=j_series+ind2;
1354        check_bounds=j1_series<MinIndex_j | j2_series>MaxIndex_j;
1355        NomTypeNc='_1_1-2';
1356    otherwise %bursts
1357        i1_series=i_series(1,:);% do not sweep the j index
1358        i2_series=i_series(1,:);
1359        j1_series=ind1*ones(1,size(i_series,2));% j index is fixed by pair choice
1360        j2_series=ind2*ones(1,size(i_series,2));
1361        check_bounds=zeros(size(i1_series));% no limitations due to min-max indices
1362end
1363
1364
1365
1366
Note: See TracBrowser for help on using the repository browser.