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

Last change on this file since 858 was 858, checked in by sommeria, 9 years ago

bugs_corrected

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