source: trunk/src/series/check_peaklock.m @ 603

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

commit the series fcts with the new form

File size: 25.8 KB
Line 
1% 'check_peaklocking': estimte peaklocking error in a civ field series TODO: UPDATE
2%------------------------------------------------------------------------
3% function ParamOut=check_peaklocking(Param)
4%
5%%%%%%%%%%% GENERAL TO ALL SERIES ACTION FCTS %%%%%%%%%%%%%%%%%%%%%%%%%%%
6%
7% This function is used in four modes by the GUI series:
8%           1) config GUI: with no input argument, the function determine the suitable GUI configuration
9%           2) interactive input: the function is used to interactively introduce input parameters, and then stops
10%           3) RUN: the function itself runs, when an appropriate input  structure Param has been introduced.
11%           4) BATCH: the function itself proceeds in BATCH mode, using an xml file 'Param' as input.
12%
13% This function is used in four modes by the GUI series:
14%           1) config GUI: with no input argument, the function determine the suitable GUI configuration
15%           2) interactive input: the function is used to interactively introduce input parameters, and then stops
16%           3) RUN: the function itself runs, when an appropriate input  structure Param has been introduced.
17%           4) BATCH: the function itself proceeds in BATCH mode, using an xml file 'Param' as input.
18%
19%OUTPUT
20% GUI_input=list of options in the GUI series.fig needed for the function
21%
22%INPUT:
23% In run mode, the input parameters are given as a Matlab structure Param copied from the GUI series.
24% In batch mode, Param is the name of the corresponding xml file containing the same information
25% In the absence of input (as activated when the current Action is selected
26% in series), the function ouput GUI_input set the activation of the needed GUI elements
27%
28% Param contains the elements:(use the menu bar command 'export/GUI config' in series to see the current structure Param)
29%    .InputTable: cell of input file names, (several lines for multiple input)
30%                      each line decomposed as {RootPath,SubDir,Rootfile,NomType,Extension}
31%    .OutputSubDir: name of the subdirectory for data outputs
32%    .OutputDirExt: directory extension for data outputs
33%    .Action: .ActionName: name of the current activated function
34%             .ActionPath:   path of the current activated function
35%    .IndexRange: set the file or frame indices on which the action must be performed
36%    .FieldTransform: .TransformName: name of the selected transform function
37%                     .TransformPath:   path  of the selected transform function
38%                     .TransformHandle: corresponding function handle
39%    .InputFields: sub structure describing the input fields withfields
40%              .FieldName: name of the field
41%              .VelType: velocity type
42%              .FieldName_1: name of the second field in case of two input series
43%              .VelType_1: velocity type of the second field in case of two input series
44%    .ProjObject: %sub structure describing a projection object (read from ancillary GUI set_object)
45%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46function  ParamOut=check_peaklocking(Param)
47
48%% set the input elements needed on the GUI series when the action is selected in the menu ActionName
49if ~exist('Param','var') % case with no input parameter
50    ParamOut={'AllowInputSort';'off';...% allow alphabetic sorting of the list of input files (options 'off'/'on', 'off' by default)
51        'WholeIndexRange';'off';...% prescribes the file index ranges from min to max (options 'off'/'on', 'off' by default)
52        'NbSlice';'on'; ...%nbre of slices ('off' by default)
53        'VelType';'two';...% menu for selecting the velocity type (options 'off'/'one'/'two',  'off' by default)
54        'FieldName';'off';...% menu for selecting the field (s) in the input file(options 'off'/'one'/'two', 'off' by default)
55        'FieldTransform'; 'off';...%can use a transform function
56        'ProjObject';'on';...%can use projection object(option 'off'/'on',
57        'Mask';'off';...%can use mask option   (option 'off'/'on', 'off' by default)
58        'OutputDirExt';'.pklock';...%set the output dir extension
59               ''};
60        return
61end
62
63%%%%%%%%%%%% STANDARD PART  %%%%%%%%%%%%
64%% select different modes,  RUN, parameter input, BATCH
65% BATCH  case: read the xml file for batch case
66if ischar(Param)
67        Param=xml2struct(Param);
68        checkrun=0;
69% RUN case: parameters introduced as the input structure Param
70else
71    hseries=guidata(Param.hseries);%handles of the GUI series
72    if isfield(Param,'Specific')&& strcmp(Param.Specific,'?')
73        checkrun=1;% will only search interactive input parameters (preparation of BATCH mode)
74    else
75        checkrun=2; % indicate the RUN option is used
76    end
77end
78ParamOut=Param; %default output
79OutputDir=[Param.OutputSubDir Param.OutputDirExt];
80
81%% root input file(s) and type
82RootPath=Param.InputTable(:,1);
83RootFile=Param.InputTable(:,3);
84SubDir=Param.InputTable(:,2);
85NomType=Param.InputTable(:,4);
86FileExt=Param.InputTable(:,5);
87[filecell,i1_series,i2_series,j1_series,j2_series]=get_file_series(Param);
88%%%%%%%%%%%%
89% The cell array filecell is the list of input file names, while
90% filecell{iview,fileindex}:
91%        iview: line in the table corresponding to a given file series
92%        fileindex: file index within  the file series,
93% i1_series(iview,ref_j,ref_i)... are the corresponding arrays of indices i1,i2,j1,j2, depending on the input line iview and the two reference indices ref_i,ref_j
94% i1_series(iview,fileindex) expresses the same indices as a 1D array in file indices
95%%%%%%%%%%%%
96NbSlice=1;%default
97if isfield(Param.IndexRange,'NbSlice')&&~isempty(Param.IndexRange.NbSlice)
98    NbSlice=Param.IndexRange.NbSlice;
99end
100nbview=1;%number of input file series (lines in InputTable)
101nbfield_j=size(i1_series{1},1); %nb of fields for the j index (bursts or volume slices)
102nbfield_i=size(i1_series{1},2); %nb of fields for the i index
103nbfield=nbfield_j*nbfield_i; %total number of fields
104nbfield_i=floor(nbfield/NbSlice);%total number of  indexes in a slice (adjusted to an integer number of slices)
105nbfield=nbfield_i*NbSlice; %total number of fields after adjustement
106
107%determine the file type on each line from the first input file
108ImageTypeOptions={'image','multimage','mmreader','video'};
109NcTypeOptions={'netcdf','civx','civdata'};
110for iview=1:nbview
111    if ~exist(filecell{iview,1}','file')
112        displ_uvmat('ERROR',['the first input file ' filecell{iview,1} ' does not exist'],checkrun)
113        return
114    end
115    [FileType{iview},FileInfo{iview},MovieObject{iview}]=get_file_type(filecell{iview,1});
116    CheckImage{iview}=~isempty(find(strcmp(FileType{iview},ImageTypeOptions)));% =1 for images
117    CheckNc{iview}=~isempty(find(strcmp(FileType{iview},NcTypeOptions)));% =1 for netcdf files
118    if ~isempty(j1_series{iview})
119        frame_index{iview}=j1_series{iview};
120    else
121        frame_index{iview}=i1_series{iview};
122    end
123end
124
125%% calibration data and timing: read the ImaDoc files
126%none
127
128%% coordinate transform or other user defined transform
129% none
130
131%%%%%%%%%%%% END STANDARD PART  %%%%%%%%%%%%
132 % EDIT FROM HERE
133
134%% check the validity of  ctinput file types
135%none
136
137%% Set field names and velocity types
138InputFields{1}=[];%default (case of images)
139if isfield(Param,'InputFields')
140    InputFields{1}=Param.InputFields;
141end
142% only one input fieldseries
143
144%% Initiate output fields
145%initiate the output structure as a copy of the first input one (reproduce fields)
146[DataOut,tild,errormsg] = read_field(filecell{1,1},FileType{1},InputFields{1},1);
147if ~isempty(errormsg)
148    displ_uvmat('ERROR',['error reading ' filecell{1,1} ': ' errormsg],checkrun)
149    return
150end
151time_1=[];
152if isfield(DataOut,'Time')
153    time_1=DataOut.Time(1);
154end
155if CheckNc{iview}
156    if isempty(strcmp('Conventions',DataOut.ListGlobalAttribute))
157        DataOut.ListGlobalAttribute=['Conventions' DataOut.ListGlobalAttribute];
158    end
159    DataOut.Conventions='uvmat';
160    DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {Param.Action}];
161    ActionKey='Action';
162    while isfield(DataOut,ActionKey)
163        ActionKey=[ActionKey '_1'];
164    end
165    DataOut.(ActionKey)=Param.Action;
166    DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {ActionKey}];
167    if isfield(DataOut,'Time')
168        DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {'Time','Time_end'}];
169    end
170end
171
172%%%%%%%%%%%%%%%% loop on field indices %%%%%%%%%%%%%%%%
173index_slice=1:nbfield;% select the file indices
174for index=index_slice
175    if checkrun
176        update_waitbar(hseries.Waitbar,index/(nbfield))
177        stopstate=get(hseries.RUN,'BusyAction');
178    else
179        stopstate='queue';
180    end
181    if isequal(stopstate,'queue')% enable STOP command
182        Data=cell(1,nbview);%initiate the set Data;
183        nbtime=0;
184        dt=[];
185        %%%%%%%%%%%%%%%% loop on views (input lines) %%%%%%%%%%%%%%%%
186        for iview=1:nbview
187            % reading input file(s)
188            [Data{iview},tild,errormsg] = read_field(filecell{iview,index},FileType{iview},InputFields{iview},frame_index{iview}(index));
189            if ~isempty(errormsg)
190                errormsg=['time_series / read_field / ' errormsg];
191                display(errormsg)
192                break
193            end
194        end
195        if isempty(errormsg)
196            Field=Data{1}; % default input field structure
197            % coordinate transform (or other user defined transform)
198            % none
199           
200            %field projection on an object
201            if Param.CheckObject
202                [Field,errormsg]=proj_field(Field,Param.ProjObject);
203                if ~isempty(errormsg)
204                    msgbox_uvmat('ERROR',['time_series / proj_field / ' errormsg])
205                    return
206                end
207            end
208           
209            % initiate the time series at the first iteration
210            if index==1
211                % stop program if the first field reading is in error
212                if ~isempty(errormsg)
213                    displ_uvmat('ERROR',['time_series / sub_field / ' errormsg],checkrun)
214                    return
215                end
216                DataOut=Field;%default
217                DataOut.NbDim=Field.NbDim+1; %add the time dimension for plots
218                nbvar=length(Field.ListVarName);
219                if nbvar==0
220                    displ_uvmat('ERROR','no input variable selected',checkrun)
221                    return
222                end
223                testsum=2*ones(1,nbvar);%initiate flag for action on each variable
224                if isfield(Field,'VarAttribute') % look for coordinate and flag variables
225                    for ivar=1:nbvar
226                        if length(Field.VarAttribute)>=ivar && isfield(Field.VarAttribute{ivar},'Role')
227                            var_role=Field.VarAttribute{ivar}.Role;%'role' of the variable
228                            if isequal(var_role,'errorflag')
229                                displ_uvmat('ERROR','do not handle error flags in time series',checkrun)
230                                return
231                            end
232                            if isequal(var_role,'warnflag')
233                                testsum(ivar)=0;  % not recorded variable
234                                eval(['DataOut=rmfield(DataOut,''' Field.ListVarName{ivar} ''');']);%remove variable
235                            end
236                            if isequal(var_role,'coord_x')| isequal(var_role,'coord_y')|...
237                                    isequal(var_role,'coord_z')|isequal(var_role,'coord')
238                                testsum(ivar)=1; %constant coordinates, record without time evolution
239                            end
240                        end
241                        % check whether the variable ivar is a dimension variable
242                        DimCell=Field.VarDimName{ivar};
243                        if ischar(DimCell)
244                            DimCell={DimCell};
245                        end
246                        if numel(DimCell)==1 && isequal(Field.ListVarName{ivar},DimCell{1})%detect dimension variables
247                            testsum(ivar)=1;
248                        end
249                    end
250                end
251                for ivar=1:nbvar
252                    if testsum(ivar)==2
253                        eval(['DataOut.' Field.ListVarName{ivar} '=[];'])
254                    end
255                end
256                DataOut.ListVarName=[{'Time'} DataOut.ListVarName];
257            end
258           
259            % add data to the current field
260            for ivar=1:length(Field.ListVarName)
261                VarName=Field.ListVarName{ivar};
262                VarVal=Field.(VarName);
263                if testsum(ivar)==2% test for recorded variable
264                    if isempty(errormsg)
265                        if isequal(Param.ProjObject.ProjMode,'inside')% take the average in the domain for 'inside' mode
266                            if isempty(VarVal)
267                                displ_uvmat('ERROR',['empty result at frame index ' num2str(i1_series{iview}(index))],checkrun)
268                                return
269                            end
270                            VarVal=mean(VarVal,1);
271                        end
272                        VarVal=shiftdim(VarVal,-1); %shift dimension
273                        DataOut.(VarName)=cat(1,DataOut.(VarName),VarVal);%concanete the current field to the time series
274                    else
275                        DataOut.(VarName)=cat(1,DataOut.(VarName),0);% put each variable to 0 in case of input reading error
276                    end
277                elseif testsum(ivar)==1% variable representing fixed coordinates
278                    VarInit=DataOut.(VarName);
279                    if isempty(errormsg) && ~isequal(VarVal,VarInit)
280                        displ_uvmat('ERROR',['time series requires constant coordinates ' VarName],checkrun)
281                        return
282                    end
283                end
284            end
285
286        end
287    end
288end
289%%%%%%% END OF LOOP WITHIN A SLICE
290
291%remove time for global attributes if exists
292Time_index=find(strcmp('Time',DataOut.ListGlobalAttribute));
293if ~isempty(Time_index)
294    DataOut.ListGlobalAttribute(Time_index)=[];
295end
296DataOut.Conventions='uvmat';
297for ivar=1:numel(DataOut.ListVarName)
298    VarName=DataOut.ListVarName{ivar};
299    eval(['DataOut.' VarName '=squeeze(DataOut.' VarName ');']) %remove singletons
300end
301
302% add time dimension
303for ivar=1:length(Field.ListVarName)
304    DimCell=Field.VarDimName(ivar);
305    if testsum(ivar)==2%variable used as time series
306        DataOut.VarDimName{ivar}=[{'Time'} DimCell];
307    elseif testsum(ivar)==1
308        DataOut.VarDimName{ivar}=DimCell;
309    end
310end
311indexremove=find(~testsum);
312if ~isempty(indexremove)
313    DataOut.ListVarName(1+indexremove)=[];
314    DataOut.VarDimName(indexremove)=[];
315    if isfield(DataOut,'Role') && ~isempty(DataOut.Role{1})%generaliser aus autres attributs
316        DataOut.Role(1+indexremove)=[];
317    end
318end
319
320%shift variable attributes
321if isfield(DataOut,'VarAttribute')
322    DataOut.VarAttribute=[{[]} DataOut.VarAttribute];
323end
324DataOut.VarDimName=[{'Time'} DataOut.VarDimName];
325DataOut.Action=Param.Action;%name of the processing programme
326test_time=diff(DataOut.Time)>0;% test that the readed time is increasing (not constant)
327if ~test_time
328    DataOut.Time=1:filecounter;
329end
330
331% display nbmissing
332if ~isequal(nbmissing,0)
333    displ_uvmat('WARNING',[num2str(nbmissing) ' files skipped: missing files or bad input, see command window display'],checkrun)
334end
335
336%name of result file
337OutputFile=fullfile_uvmat(RootPath{1},OutputDir,RootFile{1},FileExtOut,NomTypeOut,i1_series{1}(1),i1_series{1}(end),i_slice,[]);
338errormsg=struct2nc(OutputFile,DataOut); %save result file
339if isempty(errormsg)
340    display([OutputFile ' written'])
341else
342    displ_uvmat('ERROR',['error in Series/struct2nc: ' errormsg],checkrun)
343end
344
345return
346
347%%%%%%%%%%%%%%%%%%  END%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
348%evaluation of peacklocking errors
349%use splinhist: give spline coeff cc for a smooth histo (call spline4)
350%use histsmooth(x,cc): calculate the smooth histo for any value x
351%use histder(x,cc): calculate the derivative of the smooth histo
352global hfig1 hfig2 hfig3
353global nbb Uval Vval Uhist Vhist % nbb resolution of the histogram nbb=10: 10 values in unity interval
354global xval xerror yval yerror
355
356set(handles.vector_y,'Value',1)% trigger the option Uhist on the interface
357set(handles.Vhist_input,'Value',1)
358set(handles.cm_switch,'Value',0) % put the switch to 'pixel'
359
360%adjust the extremal values of the histogram in U with respect to integer
361%values
362minimU=round(min(Uval)-0.5)+0.5; %first value of the histogram with integer bins
363maximU=round(max(Uval)-0.5)+0.5;
364minim_fin=(minimU-0.5+1/(2*nbb)); % first bin valueat the beginning of an integer interval
365maxim_fin=(maximU+0.5-1/(2*nbb)); % last integer value
366nb_bin_min= round(-(minim_fin - min(Uval))*nbb); % nbre of bins added below
367nb_bin_max=round((maxim_fin -max(Uval))*nbb); %nbre of bins added above
368Uval=[minim_fin:(1/nbb):maxim_fin];
369histu_min=zeros(nb_bin_min,1);
370histu_max=zeros(nb_bin_max,1);
371Uhist=[histu_min; Uhist ;histu_max]; % column vector
372
373%adjust the extremal values of the histogram in V
374minimV=round(min(Vval-0.5)+0.5);
375maximV=round(max(Vval-0.5)+0.5);
376minim_fin=minimV-0.5+1/(2*nbb); % first bin valueat the beginning of an integer interval
377maxim_fin=maximV+0.5-1/(2*nbb); % last integer value
378nb_bin_min=round((min(Vval) - minim_fin)*nbb); % nbre of bins added below
379nb_bin_max=round((maxim_fin -max(Vval))*nbb);
380Vval=[minim_fin:(1/nbb):maxim_fin];
381histu_min=zeros(nb_bin_min,1);
382histu_max=zeros(nb_bin_max,1);
383Vhist=[histu_min; Vhist ;histu_max]; % column vector
384
385% RUN_histo_Callback(hObject, eventdata, handles)
386% %adjust the histogram to integer values:
387
388%histoU and V
389[Uhistinter,xval,xerror]=peaklock(nbb,minimU,maximU,Uhist);
390[Vhistinter,yval,yerror]=peaklock(nbb,minimV,maximV,Vhist);
391
392% selection of value ranges such that histo>=10 (enough statistics)
393Uval_ind=find(Uhist>=10);
394ind_min=min(Uval_ind);
395ind_max=max(Uval_ind);
396U_min=Uval(ind_min);% minimum allowed value
397U_max=Uval(ind_max);%maximum allowed value
398
399% selection of value ranges such that histo>=10 (enough statistics)
400Vval_ind=find(Vhist>=10);
401ind_min=min(Vval_ind);
402ind_max=max(Vval_ind);
403V_min=Vval(ind_min);% minimum allowed value
404V_max=Vval(ind_max);%maximum allowed value
405
406figure(4)% plot U histogram with smoothed one
407plot(Uval,Uhist,'b')
408grid on
409hold on
410plot(Uval,Uhistinter,'r');
411hold off
412
413figure(5)% plot V histogram with smoothed one
414plot(Vval,Vhist,'b')
415grid on
416hold on
417plot(Vval,Vhistinter,'r');
418hold off
419
420figure(6)% plot pixel error in two subplots
421hfig4=subplot(2,1,1);
422hfig5=subplot(2,1,2);
423axes(hfig4)
424plot(xval,xerror)
425axis([U_min U_max -0.4 0.4])
426xlabel('velocity u (pix)')
427ylabel('peaklocking error (pix)')
428grid on
429axes(hfig5)
430plot(yval,yerror)
431axis([V_min V_max -0.4 0.4]);
432xlabel('velocity v (pix)')
433ylabel('peaklocking error (pix)')
434grid on
435
436
437
438
439
440
441
442
443
444%'peaklock': determines peacklocking errors from velocity histograms.
445%-------------------------------------------------------
446%first smooth the input histogram 'histu' in such a way that the integral over
447%n-n+1 is preserved, then deduce the peaklocking 'error' function of the pixcel displacement 'x'.
448%
449% [histinter,x,error]=peaklock(nbb,minim,maxim,histu)
450%OUTPUT:
451%histinter: smoothed interpolated histogram
452% x: vector of displacement values.
453% error: vector of estimated errors corresponding to x
454%INPUT:
455%histu=vector representing the values of histogram  of measured velocity ;
456%minim, maxim: extremal values of the measured velocity (absica for histu)
457%nbb: number of bins inside each integer interval for the histograms
458%SUBROUTINES INCLUDED:
459%spline4.m% spline interpolation at 4th order
460%splinhist.m: give spline coeff cc for a smooth histo (call spline4)
461%histsmooth.m(x,cc): calculate the smooth histo for any value x
462%histder.m(x,cc): calculate the derivative of the smooth histo
463function [histinter,x,error]=peaklock(nbb,minim,maxim,histu)
464
465nint=maxim-minim+1
466xfin=[minim-0.5+1/(2*nbb):(1/nbb):maxim+0.5-(1/(2*nbb))];
467histo=(reshape(histu,nbb,nint));%extract values with x between integer -1/2 integer +1/2
468Integ=sum(histo)/nbb; %integral of the pdf on each integer bin
469[histinter,cc]=splinhist(Integ,minim,nbb);
470histx=reshape(histinter,nbb,nint);
471xint=[minim:1:maxim];
472x=zeros(nbb,nint);
473%determination of the displacement x(j,:)
474%j=1
475delx=histo(1,:)./histsmooth(-0.5*ones(1,nint),cc)/nbb;
476%del(1,:)=delx;
477x(1,:)=-0.5+delx-(delx.*delx/2).*histder(-0.5*ones(1,nint),cc);
478%histx(1,:)=histsmooth(x(j-1,:),cc);
479for j=2:nbb
480    delx=histo(j,:)./histsmooth(x(j-1,:),cc)/nbb;
481    %delx=delx.*(delx<3*ones(1,nint)/nbb)+3*ones(1,nint)/nbb.*~(delx <3*ones(1,nint)/nbb)
482    x(j,:)=x(j-1,:)+delx-(delx.*delx/2).*histder(x(j-1,:),cc);
483end
484%reshape
485xint=ones(nbb,1)*xint;
486x=x+xint;
487x=reshape(x,1,nbb*nint);
488error=xfin+1/(2*nbb)-x;
489
490%-------------------------------------------------------
491% --- determine the spline coefficients cc for the interpolated histogram.
492%-------------------------------------------------
493function [histsmooth,cc]= splinhist(Integ,mini,nbb)
494% provides a smooth histogramm histmooth, which remains always positive,
495% and is such that its sum over each integer bin [i-1/2 i+1/2] is equal to
496% Integ(i). The function determines histmooth as the exponential of a 4th
497% order spline function and adjust the cefficients by a Newton method to
498% fit the integral conditions Integ
499% histmooth is determined at the abscissa
500% xfin=[mini-0.5+1/(2*n):(1/n):maxi+0.5-(1/(2*n))] (maxi=mini+size(aa)-1)
501%cc(1-5,i) provides the spline coefficients
502
503% order 0
504siz=size(Integ);
505nint=siz(2);
506izero=find(Integ==0); %indices of zero elements
507inonzero=find(Integ);
508Integ(izero)=min(Integ(inonzero));
509aa=log(Integ);%initial guess for a coeff
510spli=spline4(aa,mini,nbb);  %appel à la fonction spline4
511histsmooth=exp(spli);
512
513S=(sum(reshape(histsmooth,nbb,nint)))/nbb;% integral of the fit histsmooth on ]i-1/2 i+1/2[
514epsilon=max(abs(Integ-S));
515iter=0;
516while epsilon > 0.000001 & iter<10
517ident=eye(nint);
518dSda=ones(nint);
519for j=1:nint% determination of the jacobian matrix dSda
520dhistda=spline4(ident(j,:),mini,nbb);
521expdhistda=dhistda.*histsmooth;
522dSda(j,:)=(sum(reshape(expdhistda,nbb,nint)))/nbb;
523end
524aa=aa+(Integ-S)*inv(dSda);%new estimate of coefficients aa by linear interpolation
525[spli,bb]=spline4(aa,mini,nbb);% new fit histsmooth
526histsmooth=exp(spli);
527S=(sum(reshape(histsmooth,nbb,nint)))/nbb;% integral of the fit histsmooth on ]i-1/2 i+1/2[
528epsilon=max(abs(Integ-S));
529iter=iter+1;
530end
531if iter==10, errordlg('splinhist did not converge after 10 iterations'),end
532cc(1,:)=aa;
533cc(2,:)=bb(1,:);
534cc(3,:)=bb(2,:);
535cc(4,:)=bb(3,:);
536cc(5,:)=bb(4,:);
537
538%-------------------------------------------------------
539% --- determine the 4th order spline coefficients from the function values aa.
540%-------------------------------------------------
541function [histsmooth,bb]= spline4(aa,mini,n)
542% spline interpolation at 4th order
543%aa=vector of values of a function at integer abscissa, starting at mini
544%n=number of subdivisions for the interpolated function
545% histmooth =interpolated values at absissa
546% xfin=[mini-0.5+1/(2*n):(1/n):maxi+0.5-(1/(2*n))] (maxi=mini+size(aa)-1)
547%bb=[b(i);c(i);d(i); e(i)] matrix of spline coeff
548L1=[1/2 1/4 1/8 1/16;1 1 3/4 1/2;0 2 3 3;0 0 6 12];
549L2=[-1/2 1/4 -1/8 1/16;1 -1 3/4 -1/2;0 2 -3 3;0 0 6 -12];
550M=inv(L2)*L1;
551[V,D]=eig(M);
552F=-inv(V)*inv(L2)*[1 ;0 ;0;0];
553a1rev=[1 -1/D(1,1)];
554b1rev=[F(1)/D(1,1)];
555a2rev=[1 -1/D(2,2)];
556b2rev=[F(2)/D(2,2)];
557a3=[1 -D(3,3)];
558b3=[F(3)];
559a4=[1 -D(4,4)];
560b4=[F(4)];
561
562%data
563% n=10;% résolution de la pdf: nbre de points par unite de u
564% mini=-10.0;%general mini=uint16(min(values)-1 CHOOSE maxi-mini+1 EVEN
565% maxi=9.0; % general maxi=uint16(max(values))+1
566%nint=double(maxi-mini+1); % nombre d'intervals entiers EVEN!
567siz=size(aa);
568nint=siz(2);
569maxi=mini+nint-1;
570npdf=nint*n;% nbre total d'intervals à introduire dans la pdf: hist(u,npdf)
571%simulation de pdf
572xfin=[mini-0.5+1/(2*n):(1/n):maxi+0.5-(1/(2*n))];% valeurs d'interpolation: we take n values in each integer interval
573%histolin=exp(-(xfin-1).*(xfin-1)).*(2+cos(10*(xfin-1)));% simulation d'une pdf
574%histo=log(histolin);
575%histo=sin(2*pi*xfin);
576%histextract=(reshape(histo,n,nint));
577%aa=sum(histextract)/n %integral of the pdf on each integer bin
578IP=[0 diff(aa)];
579Irev=zeros(size(aa));
580for i=1:nint
581    Irev(i)=aa(end-i+1);
582end
583IPrev=[0 diff(Irev)];
584
585%get the spline coelfficients a_d, using filter on the eigen vectors A,B,C
586Arev=filter(b1rev,a1rev,IPrev);
587Brev=filter(b2rev,a2rev,IPrev);
588C=filter(b3,a3,IP);
589D=filter(b4,a4,IP);
590A=zeros(size(Arev));
591B=zeros(size(Brev));
592for i=1:nint
593    A(i)=Arev(end-i+1);
594    B(i)=Brev(end-i+1);
595end
596%Matr=V*[A;B;C;D];
597bb=V*[A;B;C;D];
598%b=Matr(1,:);
599%c=Matr(2,:);
600%d=Matr(3,:);
601%e=Matr(4,:);
602%a=aa;
603
604%calculate the interpolation using the spline coefficients a-d
605%xextract=(reshape(xfin,n,nint));%
606chi=xfin+1/(2*n)-min(xfin)-double(int16(xfin+(1/(2*n))-min(xfin)))-0.5;% decimal part
607chi2=chi.*chi;
608chi3=chi2.*chi;
609chi4=chi3.*chi;
610avec=reshape(ones(n,1)*aa,1,n*nint);
611bvec=reshape(ones(n,1)*bb(1,:),1,n*nint);
612cvec=reshape(ones(n,1)*bb(2,:),1,n*nint);
613dvec=reshape(ones(n,1)*bb(3,:),1,n*nint);
614evec=reshape(ones(n,1)*bb(4,:),1,n*nint);
615histsmooth=avec+bvec.*chi+cvec.*chi2+dvec.*chi3+evec.*chi4;
616
617%-------------------------------------------------------
618% --- determine the interpolated histogram at points chi from the spline ceff cc.
619%-------------------------------------------------
620function histx= histsmooth(chi,cc)
621% provides the value of the interpolated histogram at values chi=x-i
622%(difference with the mnearest integer)
623% cc(5,size(chi)) is the set of spline coefficients obtained by splinhist
624chi2=chi.*chi;
625chi3=chi2.*chi;
626chi4=chi3.*chi;
627histx=exp(cc(1,:)+cc(2,:).*chi+cc(3,:).*chi2+cc(4,:).*chi3+cc(5,:).*chi4);
628
629%-------------------------------------------------------
630% --- determine the derivative p'/p of the interpolated histogram at points chi from the spline ceff cc.
631%-------------------------------------------------
632function histder= histder(chi,cc)
633% provides the logarithmique derivative p'/p of the interpolated histogram
634%at values chi=x-i
635%(difference with the nearest integer)
636% cc(5,size(chi)) is the set of spline coefficients obtained by splinhist
637chi2=chi.*chi;
638chi3=chi2.*chi;
639chi4=chi3.*chi;
640histder=cc(2,:)+2*cc(3,:).*chi+3*cc(4,:).*chi2+4*cc(5,:).*chi3;
Note: See TracBrowser for help on using the repository browser.