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

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