[Paraview] UPDATE_TIME_INDEX not part of Paraview 2.6

Sean Ziegeler seanzig at users.sourceforge.net
Wed Feb 7 13:03:40 EST 2007


Are the time steps set to 100 at some point before they get set to the 
correct number?  I've seen problems in the vtkCSCSNetCDF reader for PV 
2.4.4 (haven't tried it with other readers or versions yet) such that 
when the timesteps are set once, they can't be changed later.

For example, in vtkCSCSNetCDF, if it originally detects a single file 
with one time step and you add the wildcards to locate more files with 
additional time steps, I can't make it change it from 1 time step. 
However, if I modify the reader to notice the other files initially and 
set the correct number of time steps the first time, that works fine.

Since the above work-around has gotten me by, I haven't mentioned it 
until now.  Also, I thought it may not be an issue in 2.6.  Am I (we) 
doing something wrong, or is this a bug?

-Sean

Mike Jackson wrote:
> Yep,
>    In RequestInformation, I am setting the values; Is that the proper 
> place to set them? If I do it in RequestData then it is probably too 
> late. I have the range set from 0 to 1 in the constructor of my class.
> 
> 
> That post was by me.. I have never been able to fix the problem.
> 
> --Mike Jackson   Senior Research Engineer
> Innovative Management & Technology Services
> 
> 
> On Feb 7, 2007, at 12:13 PM, John Biddiscombe wrote:
> 
>> Mike
>>
>> make sure you also set
>>
>>    this->TimeStepRange[0] = 0;
>>    this->TimeStepRange[1] = t>0 ? this->TimeStepValues.size()-1 : 0;
>>
>> as the paraview gui will query GetTimeStepRange to find how many steps 
>> are there if this is missing, then the number of steps will be wrong.
>> Wasn't there a similar post about only seeing 100 steps recently?
>> JB
>>
>>
>>> Thanks for the code. I downloaded your OpenDX reader to get an idea 
>>> of the context. I have updated my code and things seem to work now. 
>>> The only problem that I am having is the "Time Step" slider is not 
>>> correctly updated when I open my file, it always defaults to 0->100 
>>> timesteps, regardless of what is in the file.
>>>
>>>   In my client side xml I have the following:
>>>
>>>     <Module name="PrdsDislocationReader"
>>>           root_name="PrdsDislocationReader"
>>>           output="vtkUnstructuredGrid"
>>>           module_type="Reader"
>>>           extensions=".dat .disl"
>>>           file_description="Paradis Dislocation Data">
>>>
>>>     <Source class="vtkPrdsDislocationReader"/>
>>>
>>>     <Scale property="TimeStep"
>>>            label="Time step"
>>>            trace_name="TimeStep"
>>>            keeps_timesteps="1"
>>>            help="Set the current timestep."/>
>>>     </Module>
>>>
>>> and in my server manager file I have:
>>>
>>>     <SourceProxy name="PrdsDislocationReader" 
>>> class="vtkPrdsDislocationReader">
>>>      <StringVectorProperty
>>>         name="FileName"
>>>         command="SetFileName"
>>>         number_of_elements="1">
>>>         <StringListDomain name="files"/>
>>>      </StringVectorProperty>
>>>
>>>      <IntVectorProperty
>>>         name="TimeStep"
>>>         command="SetTimeStep"
>>>         number_of_elements="1"
>>>         animateable="1"
>>>         default_values="0"
>>>         information_property="TimestepValues">
>>>
>>>        <IntRangeDomain name="range">
>>>           <RequiredProperties>
>>>              <Property name="TimeStepRangeInfo" function="Range"/>
>>>           </RequiredProperties>
>>>        </IntRangeDomain>
>>>      </IntVectorProperty>
>>>
>>>      <IntVectorProperty
>>>         name="TimeStepRangeInfo"
>>>         command="GetTimeStepRange"
>>>         information_only="1">
>>>         <SimpleIntInformationHelper/>
>>>      </IntVectorProperty>
>>>
>>>
>>>
>>>      <DoubleVectorProperty
>>>         name="TimestepValues"
>>>         information_only="1">
>>>         <TimeStepsInformationHelper/>
>>>      </DoubleVectorProperty>
>>>    <!-- End ParadisReader -->
>>>    </SourceProxy>
>>>
>>> I have tried to mimic the ExodusReader code where I can but nothing 
>>> seems to work. Does anyone have any ideas what I need to 
>>> set/over-ride or implement? Maybe something from the superclass 
>>> (vtkPolyDataAlgorithm) maybe?
>>>
>>> Thanks for any help
>>> --Mike Jackson   Senior Research Engineer
>>> Innovative Management & Technology Services
>>>
>>>
>>> On Feb 7, 2007, at 3:14 AM, John Biddiscombe wrote:
>>>
>>>> Mike, The time handling has changed a bit recently. TIME_INDEX is no 
>>>> longer supported (partly because you may have 100 time steps of data 
>>>> with step number 34 missing, if step 35 is requested, do they want 
>>>> the 35th of the data you actually have, or the 35th of the data you 
>>>> would have if the steps were all present). Instead, use a time value 
>>>> which is real and represents the actual time.(this also allows 
>>>> source which can produce a continuum of time values to be accessed)
>>>>
>>>> To get things right, try this
>>>>
>>>> //---------------------------------------------------------------------------- 
>>>>
>>>>
>>>> In RequestInformation, make sure you export the steps you can 
>>>> produce (often just 1.0, 2.0, 3.0.....etc)
>>>>
>>>>  outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),
>>>>               &this->TimeStepValues[0], this->TimeStepValues.size());
>>>>
>>>>    where you may have std::vector<double> TimeStepValues - with the 
>>>> values pushed back when the file is opened
>>>>
>>>> //---------------------------------------------------------------------------- 
>>>>
>>>>
>>>> In RequestData, check to see which time step was requested, this 
>>>> usually will be set in the GUI via SetTimeStep, but might be passed 
>>>> into the filter from downstream, if another filter is requesting 
>>>> data from a particular step
>>>>
>>>>  this->ActualTimeStep = this->TimeStep; // set by the gui, but we 
>>>> will override this if the pipeline requested another
>>>>
>>>>  if 
>>>> (outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS()))
>>>>    {
>>>>    double requestedTimeValue = 
>>>> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS())[0];
>>>>    // this is overkill, I use it to ensure that if a time value 
>>>> which is slightly out is requested
>>>>    // then the correct one is fetched
>>>>    this->ActualTimeStep = vtkstd::find_if(
>>>>      this->TimeStepValues.begin(),
>>>>      this->TimeStepValues.end(),
>>>>      vtkstd::bind2nd( WithinTolerance( ), requestedTimeValue ))
>>>>      - this->TimeStepValues.begin();
>>>>    this->ActualTimeStep = this->ActualTimeStep + 
>>>> this->TimeStepRange[0];
>>>>    doOutput->GetInformation()->Set(vtkDataObject::DATA_TIME_STEPS(), 
>>>> &requestedTimeValue, 1);
>>>>    CSCSOutputMacro(<<"Got a timestep request from downstream t= " << 
>>>> requestedTimeValue << " Step : " << this->ActualTimeStep);
>>>>    }
>>>>  else
>>>>    {
>>>>    double timevalue[1];
>>>>    timevalue[0] = 
>>>> this->TimeStepValues[this->ActualTimeStep-this->TimeStepRange[0]];
>>>>    CSCSOutputMacro(<<"Using manually set t= " << timevalue[0] << " 
>>>> Step : " << this->ActualTimeStep);
>>>>    doOutput->GetInformation()->Set(vtkDataObject::DATA_TIME_STEPS(), 
>>>> &timevalue[0], 1);
>>>>    }
>>>>
>>>> //---------------------------------------------------------------------------- 
>>>>
>>>> note that
>>>>
>>>> class WithinTolerance: public std::binary_function<double, double, 
>>>> bool>
>>>> {
>>>> public:
>>>>    result_type operator()(first_argument_type a, 
>>>> second_argument_type b) const
>>>>    {
>>>>      bool result = (fabs(a-b)<=(a*1E-6));
>>>>      return (result_type)result;
>>>>    }
>>>> };
>>>>
>>>> if you want to see a filter that uses the above code and works well, 
>>>> look here
>>>> https://svn.cscs.ch/vtkContrib/trunk/vtkCSCS/vtkOpenDX/
>>>> for my openDX reader which is the most simple one I've made that 
>>>> uses the above code.
>>>>
>>>> JB
>>>>> I was trying to recompile my ParaView modules and I get the error:
>>>>>
>>>>> /Users/mjackson/Task_4/Workspace/PVDislocation/vtkPrdsDislocationReader.cpp:80: 
>>>>> error: 'UPDATE_TIME_INDEX' is not a member of 
>>>>> 'vtkStreamingDemandDrivenPipeline'
>>>>>
>>>>> The offending line is:
>>>>>   
>>>>> outInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_INDEX(), 
>>>>> this->GetTimeStep() );
>>>>>
>>>>> This used to work in ParaView 2.4.4 and I vaguely remember that 
>>>>> time support was being changed in PV 2.6. Could someone explain 
>>>>> what I should be using now?
>>>>>
>>>>> I have written a custom ParaView reader module for our file format. 
>>>>> I need Time Support and I was successfully using this code in PV 2.4.
>>>>>
>>>>> Thanks
>>>>> --Mike Jackson   Senior Research Engineer
>>>>> Innovative Management & Technology Services
>>>>>
>>>>>
>>>>>
>>>>> _______________________________________________
>>>>> ParaView mailing list
>>>>> ParaView at paraview.org
>>>>> http://www.paraview.org/mailman/listinfo/paraview
>>>>
>>>>
>>>> --John Biddiscombe,                            email:biddisco @ cscs.ch
>>>> http://www.cscs.ch/about/BJohn.php
>>>> CSCS, Swiss National Supercomputing Centre  | Tel:  +41 (91) 610.82.07
>>>> Via Cantonale, 6928 Manno, Switzerland      | Fax:  +41 (91) 610.82.82
>>>>
>>>>
>>>
>>>
>>> _______________________________________________
>>> ParaView mailing list
>>> ParaView at paraview.org
>>> http://www.paraview.org/mailman/listinfo/paraview
>>
>>
>> --John Biddiscombe,                            email:biddisco @ cscs.ch
>> http://www.cscs.ch/about/BJohn.php
>> CSCS, Swiss National Supercomputing Centre  | Tel:  +41 (91) 610.82.07
>> Via Cantonale, 6928 Manno, Switzerland      | Fax:  +41 (91) 610.82.82
>>
>>
> 
> 
> _______________________________________________
> ParaView mailing list
> ParaView at paraview.org
> http://www.paraview.org/mailman/listinfo/paraview
> 


More information about the ParaView mailing list