source: trunk/src/series/time_series.m @ 1085

Last change on this file since 1085 was 1085, checked in by sommeria, 4 years ago

aver_spectral added andsmall bug repairs

File size: 26.0 KB
RevLine 
[572]1%'time_series': extract a time series after projection on an object (points , line..)
[457]2% this function can be used as a template for applying a global operation on a series of input fields
[169]3%------------------------------------------------------------------------
[374]4% function GUI_input=time_series(Param)
[169]5%
[457]6%%%%%%%%%%% GENERAL TO ALL SERIES ACTION FCTS %%%%%%%%%%%%%%%%%%%%%%%%%%%
7%
[169]8%OUTPUT
[596]9% ParamOut: sets options in the GUI series.fig needed for the function
[169]10%
11%INPUT:
[457]12% In run mode, the input parameters are given as a Matlab structure Param copied from the GUI series.
13% In batch mode, Param is the name of the corresponding xml file containing the same information
[596]14% when Param.Action.RUN=0 (as activated when the current Action is selected
15% in series), the function ouput paramOut set the activation of the needed GUI elements
[169]16%
[596]17% Param contains the elements:(use the menu bar command 'export/GUI config' in series to
18% see the current structure Param)
[457]19%    .InputTable: cell of input file names, (several lines for multiple input)
20%                      each line decomposed as {RootPath,SubDir,Rootfile,NomType,Extension}
21%    .OutputSubDir: name of the subdirectory for data outputs
[474]22%    .OutputDirExt: directory extension for data outputs
[457]23%    .Action: .ActionName: name of the current activated function
24%             .ActionPath:   path of the current activated function
[596]25%             .ActionExt: fct extension ('.m', Matlab fct, '.sh', compiled   Matlab fct
26%             .RUN =0 for GUI input, =1 for function activation
27%             .RunMode='local','background', 'cluster': type of function  use
28%             
[457]29%    .IndexRange: set the file or frame indices on which the action must be performed
30%    .FieldTransform: .TransformName: name of the selected transform function
31%                     .TransformPath:   path  of the selected transform function
32%    .InputFields: sub structure describing the input fields withfields
[596]33%              .FieldName: name(s) of the field
[457]34%              .VelType: velocity type
35%              .FieldName_1: name of the second field in case of two input series
36%              .VelType_1: velocity type of the second field in case of two input series
[596]37%              .Coord_y: name of y coordinate variable
38%              .Coord_x: name of x coordinate variable
[457]39%    .ProjObject: %sub structure describing a projection object (read from ancillary GUI set_object)
40%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
[596]41
[810]42%=======================================================================
[1071]43% Copyright 2008-2020, LEGI UMR 5519 / CNRS UGA G-INP, Grenoble, France
[810]44%   http://www.legi.grenoble-inp.fr
45%   Joel.Sommeria - Joel.Sommeria (A) legi.cnrs.fr
46%
47%     This file is part of the toolbox UVMAT.
48%
49%     UVMAT is free software; you can redistribute it and/or modify
50%     it under the terms of the GNU General Public License as published
51%     by the Free Software Foundation; either version 2 of the license,
52%     or (at your option) any later version.
53%
54%     UVMAT is distributed in the hope that it will be useful,
55%     but WITHOUT ANY WARRANTY; without even the implied warranty of
56%     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
57%     GNU General Public License (see LICENSE.txt) for more details.
58%=======================================================================
59
[457]60function ParamOut=time_series(Param)
[27]61
[606]62%% set the input elements needed on the GUI series when the action is selected in the menu ActionName or InputTable refreshed
[875]63if isstruct(Param) && isequal(Param.Action.RUN,0)% function activated from the GUI series but not RUN
[605]64    ParamOut.AllowInputSort='off';% allow alphabetic sorting of the list of input file SubDir (options 'off'/'on', 'off' by default)
65    ParamOut.WholeIndexRange='off';% prescribes the file index ranges from min to max (options 'off'/'on', 'off' by default)
66    ParamOut.NbSlice='on'; %nbre of slices ('off' by default)
67    ParamOut.VelType='two';% menu for selecting the velocity type (options 'off'/'one'/'two',  'off' by default)
68    ParamOut.FieldName='two';% menu for selecting the field (s) in the input file(options 'off'/'one'/'two', 'off' by default)
[1001]69    ParamOut.FieldTransform = {'phys','phys_polar'};%can use a transform function, proposed list (needed for compilation)
70    ParamOut.TransformPath=fullfile(fileparts(which('uvmat')),'transform_field');% path to transform functions (needed for compilation only)   
71    if 1==2 % loop used to enforce compilation of transform fct, never entered.
72        phys
73        phys_polar
74    end
[605]75    ParamOut.ProjObject='on';%can use projection object(option 'off'/'on',
76    ParamOut.Mask='off';%can use mask option   (option 'off'/'on', 'off' by default)
77    ParamOut.OutputDirExt='.tseries';%set the output dir extension
78    ParamOut.OutputFileMode='NbSlice';% '=NbInput': 1 output file per input file index, '=NbInput_i': 1 file per input file index i, '=NbSlice': 1 file per slice
[871]79    % check for selection of a projection object
80    hseries=findobj(allchild(0),'Tag','series');% handles of the GUI series
[875]81    hhseries=guidata(hseries);
[871]82    if  ~isfield(Param,'ProjObject')
[1085]83        msgbox_uvmat('WARNING','you may need to introduce a projection object for the time_series');
[871]84    end
[875]85   
[871]86    % introduce bin size for histograms
[872]87    if isfield(Param,'CheckObject') &&Param.CheckObject
[871]88        SeriesData=get(hseries,'UserData');
[872]89        if ismember(SeriesData.ProjObject.ProjMode,{'inside','outside'})
[875]90            answer=msgbox_uvmat('INPUT_TXT','set bin size for histograms (or keep ''auto'' by default)?','auto');
[880]91            ParamOut.ActionInput.VarMesh=str2num(answer);
[871]92        end
93    end
[875]94   
[874]95    % test for subtraction
[875]96    if size(Param.InputTable,1)>2
97        msgbox_uvmat('WARNING','''time_series'' uses only one or two input lines (two for substraction). To concatene fields first use ''merge_proj''');
98    end
99    if size(Param.InputTable,1)>=2
100        answer=msgbox_uvmat('INPUT_Y-N','substract the two input file series?');
101        if strcmp(answer,'Yes')
102            if isempty(Param.FieldTransform.TransformName)
103                set(hhseries.TransformName,'value',2) %select sub_field
104            end
105        else
106            set(hhseries.InputTable,'Data',Param.InputTable(1,:))
[874]107        end
[875]108    end
109   
[871]110    % check the existence of the first and last file in the series
[875]111    first_j=[];
[716]112    if isfield(Param.IndexRange,'first_j'); first_j=Param.IndexRange.first_j; end
113    last_j=[];
114    if isfield(Param.IndexRange,'last_j'); last_j=Param.IndexRange.last_j; end
115    PairString='';
116    if isfield(Param.IndexRange,'PairString'); PairString=Param.IndexRange.PairString; end
117    [i1,i2,j1,j2] = get_file_index(Param.IndexRange.first_i,first_j,PairString);
118    FirstFileName=fullfile_uvmat(Param.InputTable{1,1},Param.InputTable{1,2},Param.InputTable{1,3},...
119        Param.InputTable{1,5},Param.InputTable{1,4},i1,i2,j1,j2);
120    if ~exist(FirstFileName,'file')
121        msgbox_uvmat('WARNING',['the first input file ' FirstFileName ' does not exist'])
[871]122    else
123        [i1,i2,j1,j2] = get_file_index(Param.IndexRange.last_i,last_j,PairString);
124        LastFileName=fullfile_uvmat(Param.InputTable{1,1},Param.InputTable{1,2},Param.InputTable{1,3},...
125            Param.InputTable{1,5},Param.InputTable{1,4},i1,i2,j1,j2);
126        if ~exist(LastFileName,'file')
127            msgbox_uvmat('WARNING',['the last input file ' LastFileName ' does not exist'])
128        end
[605]129    end
[599]130    return
[27]131end
132
[526]133%%%%%%%%%%%% STANDARD PART  %%%%%%%%%%%%
[624]134ParamOut=[]; %default output
[592]135%% read input parameters from an xml file if input is a file name (batch mode)
136checkrun=1;
[457]137if ischar(Param)
[592]138    Param=xml2struct(Param);% read Param as input file (batch case)
139    checkrun=0;
[374]140end
[624]141hseries=findobj(allchild(0),'Tag','series');
142RUNHandle=findobj(hseries,'Tag','RUN');%handle of RUN button in GUI series
143WaitbarHandle=findobj(hseries,'Tag','Waitbar');%handle of waitbar in GUI series
[592]144
[624]145%% define the directory for result file (with path=RootPath{1})
[474]146OutputDir=[Param.OutputSubDir Param.OutputDirExt];
[457]147
[609]148%% root input file(s) name, type and index series
[457]149RootPath=Param.InputTable(:,1);
150RootFile=Param.InputTable(:,3);
151SubDir=Param.InputTable(:,2);
152NomType=Param.InputTable(:,4);
153FileExt=Param.InputTable(:,5);
[960]154
[751]155hdisp=disp_uvmat('WAITING...','checking the file series',checkrun);
[960]156% gives the series of input file names and indices set by the input parameters:
[374]157[filecell,i1_series,i2_series,j1_series,j2_series]=get_file_series(Param);
[526]158% filecell{iview,fileindex}:
[457]159%        iview: line in the table corresponding to a given file series
[960]160%        fileindex: file index with i and j reshaped as a 1D array
161% 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
[457]162% i1_series(iview,fileindex) expresses the same indices as a 1D array in file indices
[960]163if ~isempty(hdisp),delete(hdisp),end;%end the waiting display
[374]164
[960]165NbView=numel(i1_series);%number of input file series (lines in InputTable)
166NbField_j=size(i1_series{1},1); %nb of fields for the j index (bursts or volume slices)
167NbField_i=size(i1_series{1},2); %nb of fields for the i index
168NbField=NbField_j*NbField_i; %total number of fields
169
[609]170%% determine the file type on each line from the first input file
[984]171ImageTypeOptions={'image','multimage','mmreader','video','cine_phantom'};
[457]172NcTypeOptions={'netcdf','civx','civdata'};
[970]173FileType=cell(1,NbView);
174FileInfo=cell(1,NbView);
175MovieObject=cell(1,NbView);
176CheckImage=cell(1,NbView);
177CheckNc=cell(1,NbView);
178frame_index=cell(1,NbView);
[874]179
[970]180for iview=1:NbView
[457]181    if ~exist(filecell{iview,1}','file')
[668]182        disp_uvmat('ERROR',['the first input file ' filecell{iview,1} ' does not exist'],checkrun)
[27]183        return
184    end
[784]185    [FileInfo{iview},MovieObject{iview}]=get_file_info(filecell{iview,1});
[783]186    FileType{iview}=FileInfo{iview}.FileType;
[768]187    if strcmp(FileType{iview},'civdata')||strcmp(FileType{iview},'civx')
188        if ~isfield(Param.InputFields,'VelType')
189            FileType{iview}='netcdf';% civ data read as usual netcdf files
190        end
191    end
[457]192    CheckImage{iview}=~isempty(find(strcmp(FileType{iview},ImageTypeOptions)));% =1 for images
193    CheckNc{iview}=~isempty(find(strcmp(FileType{iview},NcTypeOptions)));% =1 for netcdf files
[768]194    if isempty(j1_series{iview})
195        frame_index{iview}=i1_series{iview};
196    else
[457]197        frame_index{iview}=j1_series{iview};
198    end
[27]199end
200
[457]201%% calibration data and timing: read the ImaDoc files
[470]202[XmlData,NbSlice_calib,time,errormsg]=read_multimadoc(RootPath,SubDir,RootFile,FileExt,i1_series,i2_series,j1_series,j2_series);
[768]203if ~isempty(errormsg)
204    disp_uvmat('ERROR',['error in reading xmlfile: ' errormsg],checkrun)
205    return
206end
[470]207if size(time,1)>1
[457]208    diff_time=max(max(diff(time)));
209    if diff_time>0
[668]210        disp_uvmat('WARNING',['times of series differ by (max) ' num2str(diff_time)],checkrun)
[599]211    end
[531]212    time=time(1,:);% choose the time data from the first sequence
[457]213end
[27]214
[457]215%% coordinate transform or other user defined transform
[526]216transform_fct=[];%default
[478]217if isfield(Param,'FieldTransform')&&~isempty(Param.FieldTransform.TransformName)
[1001]218    %addpath(Param.FieldTransform.TransformPath)
[474]219    transform_fct=str2func(Param.FieldTransform.TransformName);
[1001]220    %rmpath(Param.FieldTransform.TransformPath)
[883]221    if isfield(Param,'TransformInput')
222        XmlData{1}.TransformInput=Param.TransformInput;
223    end
[457]224end
[474]225
[457]226%%%%%%%%%%%% END STANDARD PART  %%%%%%%%%%%%
[599]227% EDIT FROM HERE
[27]228
[656]229%% check the validity of  the input file types
230if ~CheckImage{1}&&~CheckNc{1}
[668]231    disp_uvmat('ERROR',['invalid file type input ' FileType{1}],checkrun)
[457]232    return
[27]233end
[970]234if NbView==2 && ~isequal(CheckImage{1},CheckImage{2})
[668]235    disp_uvmat('ERROR','input must be two image series or two netcdf file series',checkrun)
[259]236    return
237end
[635]238
239%% settings for the output file
[656]240FileExtOut='.nc';% write result as .nc files for netcdf inputs
[609]241NomTypeOut=nomtype2pair(NomType{1});% determine the index nomenclature type for the output file
[635]242first_i=i1_series{1}(1);
243last_i=i1_series{1}(end);
244if isempty(j1_series{1})% if there is no second index j
245    first_j=1;last_j=1;
246else
247    first_j=j1_series{1}(1);
248    last_j=j1_series{1}(end);
249end
[457]250
251%% Set field names and velocity types
252InputFields{1}=[];%default (case of images)
[970]253if NbView==2
[867]254    InputFields{2}=[];%default (case of images)
255end
[457]256if isfield(Param,'InputFields')
257    InputFields{1}=Param.InputFields;
[970]258    if NbView==2
[867]259        InputFields{2}=Param.InputFields;%default
[457]260        if isfield(Param.InputFields,'FieldName_1')
261            InputFields{2}.FieldName=Param.InputFields.FieldName_1;
262            if isfield(Param.InputFields,'VelType_1')
263                InputFields{2}.VelType=Param.InputFields.VelType_1;
264            end
265        end
266    end
267end
[27]268
[457]269%% Initiate output fields
270%initiate the output structure as a copy of the first input one (reproduce fields)
[464]271[DataOut,tild,errormsg] = read_field(filecell{1,1},FileType{1},InputFields{1},1);
[457]272if ~isempty(errormsg)
[668]273    disp_uvmat('ERROR',['error reading ' filecell{1,1} ': ' errormsg],checkrun)
[457]274    return
275end
276time_1=[];
277if isfield(DataOut,'Time')
278    time_1=DataOut.Time(1);
279end
280if CheckNc{iview}
281    if isempty(strcmp('Conventions',DataOut.ListGlobalAttribute))
282        DataOut.ListGlobalAttribute=['Conventions' DataOut.ListGlobalAttribute];
283    end
284    DataOut.Conventions='uvmat';
285    DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {Param.Action}];
286    ActionKey='Action';
287    while isfield(DataOut,ActionKey)
288        ActionKey=[ActionKey '_1'];
289    end
290    DataOut.(ActionKey)=Param.Action;
291    DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {ActionKey}];
292    if isfield(DataOut,'Time')
293        DataOut.ListGlobalAttribute=[DataOut.ListGlobalAttribute {'Time','Time_end'}];
294    end
295end
296
[871]297nbfile=0;% not used , to check
[599]298nbmissing=0;
[880]299VarMesh=[];
[871]300checkhisto=0;
[1069]301checkline=0;
302if isfield(Param,'ProjObject')
303    if ismember(Param.ProjObject.ProjMode,{'inside','outside'})
[871]304    checkhisto=1;
[880]305    if isfield(Param,'ActionInput') && isfield(Param.ActionInput,'VarMesh')%case of histograms
306        VarMesh=Param.ActionInput.VarMesh;
[884]307    else
308        VarMesh=[];
309        disp_uvmat('WARNING','automatic bin size for histograms, select time_series again to set the value',checkrun)
[871]310    end
[1069]311    elseif ismember(Param.ProjObject.Type,{'line'})
312        checkline=1;
313    end
[871]314end
[875]315
[599]316%%%%%%%%%%%%%%%% loop on field indices %%%%%%%%%%%%%%%%
[970]317for index=1:NbField
318    update_waitbar(WaitbarHandle,index/NbField)
[635]319    if ~isempty(RUNHandle) && ~strcmp(get(RUNHandle,'BusyAction'),'queue')
[624]320        disp('program stopped by user')
321        break % leave the loop if stop is ordered
[599]322    end
[970]323    Data=cell(1,NbView);%initiate the set Data;
[599]324    nbtime=0;
325    dt=[];
326    %%%%%%%%%%%%%%%% loop on views (input lines) %%%%%%%%%%%%%%%%
[970]327    for iview=1:NbView
[599]328        % reading input file(s)
329        [Data{iview},tild,errormsg] = read_field(filecell{iview,index},FileType{iview},InputFields{iview},frame_index{iview}(index));
330        if ~isempty(errormsg)
331            errormsg=['time_series / read_field / ' errormsg];
[880]332            disp(errormsg)
333            break% exit the loop on iview
[442]334        end
[599]335        if ~isempty(NbSlice_calib)
336            Data{iview}.ZIndex=mod(i1_series{iview}(index)-1,NbSlice_calib{iview})+1;%Zindex for phys transform
337        end
338    end
[880]339    if ~isempty(errormsg)
340        continue %in case of input error skip the current input field, go to next index
341    end
342    Field=Data{1}; % default input field structure
343    % coordinate transform (or other user defined transform)
344    if ~isempty(transform_fct)
345        switch nargin(transform_fct)
346            case 4
347                if length(Data)==2
348                    Field=transform_fct(Data{1},XmlData{1},Data{2},XmlData{2});
349                else
[875]350                    Field=transform_fct(Data{1},XmlData{1});
[880]351                end
352            case 3
353                if length(Data)==2
354                    Field=transform_fct(Data{1},XmlData{1},Data{2});
355                else
356                    Field=transform_fct(Data{1},XmlData{1});
357                end
358            case 2
359                Field=transform_fct(Data{1},XmlData{1});
360            case 1
361                Field=transform_fct(Data{1});
[599]362        end
[880]363    end
364   
365    %field projection on an object
366    if Param.CheckObject
367        % calculate tps coefficients if needed
368        if isfield(Param.ProjObject,'ProjMode')&& strcmp(Param.ProjObject.ProjMode,'interp_tps')
369            Field=tps_coeff_field(Field,check_proj_tps);
[599]370        end
[880]371        [Field,errormsg]=proj_field(Field,Param.ProjObject,VarMesh);
372        if ~isempty(errormsg)
373            disp_uvmat('ERROR',['time_series / proj_field / ' errormsg],checkrun)
374            return
375        end
376    end
377    %         nbfile=nbfile+1;
378   
379    % initiate the time series at the first iteration
380    if index==1
381        % stop program if the first field reading is in error
382        if ~isempty(errormsg)
383            disp_uvmat('ERROR',['time_series / sub_field / ' errormsg],checkrun)
384            return
385        end
386        if ~isfield(Field,'ListVarName')
387            disp_uvmat('ERROR','no variable in the projected field',checkrun)
388            return
389        end
390        DataOut=Field;%output reproduced the first projected field by default
391        nbvar=length(Field.ListVarName);
392        if nbvar==0
393            disp_uvmat('ERROR','no input variable selected',checkrun)
394            return
395        end
396        if checkhisto%case of histograms
397            testsum=zeros(1,nbvar);%initiate flag for action on each variable
398            for ivar=1:numel(Field.ListVarName)% list of variable names before projection (histogram)
399                VarName=Field.ListVarName{ivar};
400                if isfield(Data{1},VarName)
401                    DataOut.ListVarName=[DataOut.ListVarName {[VarName 'Histo']}];
402                    DataOut.VarDimName=[DataOut.VarDimName {{'Time',VarName}}];
403                    StatName=pdf2stat;% get the names of statistical quantities to calcuilate at each time
404                    for istat=1:numel(StatName)
405                        DataOut.ListVarName=[DataOut.ListVarName {[VarName StatName{istat}]}];
406                        DataOut.VarDimName=[DataOut.VarDimName {'Time'}];
[871]407                    end
[880]408                    testsum(ivar)=1;
409                    DataOut.(VarName)=Field.(VarName);
[970]410                    DataOut.([VarName 'Histo'])=zeros([NbField numel(DataOut.(VarName))]);
[880]411                    VarMesh=Field.(VarName)(2)-Field.(VarName)(1);
[871]412                end
[880]413            end
414            disp(['mesh for histogram = ' num2str(VarMesh)])
415        else
416            testsum=2*ones(1,nbvar);%initiate flag for action on each variable
417            if isfield(Field,'VarAttribute') % look for coordinate and flag variables
418                for ivar=1:nbvar
419                    if length(Field.VarAttribute)>=ivar && isfield(Field.VarAttribute{ivar},'Role')
420                        var_role=Field.VarAttribute{ivar}.Role;%'role' of the variable
421                        if isequal(var_role,'errorflag')
422                            disp_uvmat('ERROR','do not handle error flags in time series',checkrun)
423                            return
[27]424                        end
[880]425                        if isequal(var_role,'warnflag')
426                            testsum(ivar)=0;  % not recorded variable
427                            eval(['DataOut=rmfield(DataOut,''' Field.ListVarName{ivar} ''');']);%remove variable
[27]428                        end
[1069]429                        if strcmp(var_role,'coord_x')||strcmp(var_role,'coord_z')||strcmp(var_role,'coord')
[880]430                            testsum(ivar)=1; %constant coordinates, record without time evolution
[474]431                        end
[1069]432                        if strcmp(var_role,'coord_y')&& ~checkline
433                             testsum(ivar)=1;
434                        end
[27]435                    end
[880]436                    % check whether the variable ivar is a dimension variable
437                    DimCell=Field.VarDimName{ivar};
438                    if ischar(DimCell)
439                        DimCell={DimCell};
[474]440                    end
[880]441                    if numel(DimCell)==1 && isequal(Field.ListVarName{ivar},DimCell{1})%detect dimension variables
442                        testsum(ivar)=1;
443                    end
[474]444                end
445            end
[880]446            for ivar=1:nbvar
447                if testsum(ivar)==2
448                    VarName=Field.ListVarName{ivar};
449                    siz=size(Field.(VarName));
[970]450                    DataOut.(VarName)=zeros([NbField siz]);
[330]451                end
[871]452            end
[880]453        end
454        DataOut.ListVarName=[{'Time'} DataOut.ListVarName];
455    end
456    % end initialisation for index==1
457   
458    % append data from the current field
459
460    if checkhisto % case of histogram (projection mode=inside or outside)
461        for ivar=1:length(Field.ListVarName)
462            VarName=Field.ListVarName{ivar};
463            if isfield(Data{1},VarName)
464                MaxValue=max(DataOut.(VarName));% current max of histogram absissa
465                MinValue=min(DataOut.(VarName));% current min of histogram absissa
466                MaxIndex=round(MaxValue/VarMesh);
467                MinIndex=round(MinValue/VarMesh);
468                MaxIndex_new=round(max(Field.(VarName)/VarMesh));% max of the current field
469                MinIndex_new=round(min(Field.(VarName)/VarMesh));
470                if MaxIndex_new>MaxIndex% the variable max for the current field exceeds the previous one
471                    DataOut.(VarName)=[DataOut.(VarName) VarMesh*(MaxIndex+1:MaxIndex_new)];% append the new variable values
[970]472                    DataOut.([VarName 'Histo'])=[DataOut.([VarName 'Histo']) zeros(NbField,MaxIndex_new-MaxIndex)]; % append the new histo values
[599]473                end
[880]474                if MinIndex_new <= MinIndex-1
475                    DataOut.(VarName)=[VarMesh*(MinIndex_new:MinIndex-1) DataOut.(VarName)];% insert the new variable values
[970]476                    DataOut.([VarName 'Histo'])=[zeros(NbField,MinIndex-MinIndex_new) DataOut.([VarName 'Histo'])];% insert the new histo values
[880]477                    ind_start=1;
478                else
479                    ind_start=MinIndex_new-MinIndex+1;
480                end
481                DataOut.([VarName 'Histo'])(index,ind_start:ind_start+MaxIndex_new-MinIndex_new)=...
482                    DataOut.([VarName 'Histo'])(index,ind_start:ind_start+MaxIndex_new-MinIndex_new)+Field.([VarName 'Histo']);
[881]483                VarVal=pdf2stat((Field.(VarName))',(Field.([VarName 'Histo']))');% max of the current field
[880]484                for istat=1:numel(VarVal)
485                    DataOut.([VarName StatName{istat}])(index)=VarVal(istat);
486                end
[27]487            end
[599]488        end
[880]489    else % not histogram
490        for ivar=1:length(Field.ListVarName)
491            VarName=Field.ListVarName{ivar};
492            VarVal=Field.(VarName);
493            if testsum(ivar)==2% test for recorded variable
494                if isempty(errormsg)
495                    VarVal=shiftdim(VarVal,-1); %shift dimension
496                    DataOut.(VarName)(index,:,:)=VarVal;%concanete the current field to the time series
497                end
498            elseif testsum(ivar)==1% variable representing fixed coordinates
499                VarInit=DataOut.(VarName);
500                if isempty(errormsg) && ~isequal(VarVal,VarInit)
501                    disp_uvmat('ERROR',['time series requires constant coordinates ' VarName ': use projection mode interp'],checkrun)
502                    return
503                end
[474]504            end
[27]505        end
[880]506    end
507   
508    % record the time:
509    if isempty(time)% time not set by xml filer(s)
510        if isfield(Data{1},'Time')
511            DataOut.Time(index,1)=Field.Time;
512        else
513            DataOut.Time(index,1)=index;%default
[27]514        end
[880]515    else % time from ImaDoc prevails  TODO: correct
516        DataOut.Time(index,1)=time(index);%
[330]517    end
[881]518    index
[880]519    % record the number of missing input fields
520    if ~isempty(errormsg)
521        nbmissing=nbmissing+1;
522        display(['index=' num2str(index) ':' errormsg])
523    end
[599]524end
[880]525
[599]526%%%%%%% END OF LOOP WITHIN A SLICE
527
528%remove time for global attributes if exists
529Time_index=find(strcmp('Time',DataOut.ListGlobalAttribute));
530if ~isempty(Time_index)
531    DataOut.ListGlobalAttribute(Time_index)=[];
532end
533DataOut.Conventions='uvmat';
534for ivar=1:numel(DataOut.ListVarName)
535    VarName=DataOut.ListVarName{ivar};
536    eval(['DataOut.' VarName '=squeeze(DataOut.' VarName ');']) %remove singletons
537end
538
539% add time dimension
540for ivar=1:length(Field.ListVarName)
[751]541    DimCell=Field.VarDimName{ivar};
542    if ischar(DimCell),DimCell={DimCell};end
543    if testsum(ivar)==2% variable for which time series is calculated
[599]544        DataOut.VarDimName{ivar}=[{'Time'} DimCell];
[751]545    elseif testsum(ivar)==1 % variable represneting a fixed coordinate
[599]546        DataOut.VarDimName{ivar}=DimCell;
[330]547    end
[599]548end
549indexremove=find(~testsum);
550if ~isempty(indexremove)
551    DataOut.ListVarName(1+indexremove)=[];
552    DataOut.VarDimName(indexremove)=[];
553    if isfield(DataOut,'Role') && ~isempty(DataOut.Role{1})%generaliser aus autres attributs
554        DataOut.Role(1+indexremove)=[];
[330]555    end
[27]556end
[107]557
[599]558%shift variable attributes
559if isfield(DataOut,'VarAttribute')
560    DataOut.VarAttribute=[{[]} DataOut.VarAttribute];
561end
562DataOut.VarDimName=[{'Time'} DataOut.VarDimName];
563DataOut.Action=Param.Action;%name of the processing programme
564test_time=diff(DataOut.Time)>0;% test that the readed time is increasing (not constant)
565if ~test_time
[970]566    DataOut.Time=1:NbField;
[599]567end
568
569if ~isequal(nbmissing,0)
[668]570    disp_uvmat('WARNING',[num2str(nbmissing) ' files skipped: missing files or bad input, see command window display'],checkrun)
[599]571end
572
[609]573%% name of result file
574OutputFile=fullfile_uvmat(RootPath{1},OutputDir,RootFile{1},FileExtOut,NomTypeOut,first_i,last_i,first_j,last_j);
[599]575errormsg=struct2nc(OutputFile,DataOut); %save result file
576if isempty(errormsg)
577    display([OutputFile ' written'])
578else
[668]579    disp_uvmat('ERROR',['error in Series/struct2nc: ' errormsg],checkrun)
[599]580end
581
[880]582%% plot the time series for  (the last one in case of multislices)
583if checkrun %&& isfield(Param,'ProjObject') && strcmp(Param.ProjObject.Type,'points')
584    %% open the result file with uvmat (in RUN mode)
585    uvmat(OutputFile)% open the last result file with uvmat
[474]586end
[880]587%     figure
588%     haxes=axes;
589%     plot_field(DataOut,haxes)
590%     
591%     %% display the result file using the GUI get_field
592%     hget_field=findobj(allchild(0),'name','get_field');
593%     if ~isempty(hget_field)
594%         delete(hget_field)
595%     end
596%     get_field(OutputFile,DataOut)
597% end
[27]598
599
Note: See TracBrowser for help on using the repository browser.