[vtkusers] Time Interpolation

kenichiro yoshimi rccm.kyoshimi at gmail.com
Tue May 15 20:22:22 EDT 2018


Hi Andrew,

I use the second vtkExtractTimeSteps to extract the hoping time step
from the output of vtkTemporalInterpolator
which contains multiple time steps like [0, 1(index of specfiedTime),
...] in my code. Especially, specified time step index is needed to be
"1" in extractTime2 filter:
  extractTime2.AddTimeStepIndex(1)
because the index of hoping time is made indicate "1" in
vtkTemporalInterpolator.

Best Regards

2018-05-15 23:37 GMT+09:00 Slaughter, Andrew E <andrew.slaughter at inl.gov>:
> This is great, thanks for taking time to get this working. I have one
> question, what is the purpose of the second vtkExtractTimeSteps?
>
> On Tue, May 15, 2018 at 2:13 AM, kenichiro yoshimi <rccm.kyoshimi at gmail.com>
> wrote:
>>
>> Hello Andrew,
>>
>> I found a dirty way using vtkTemporalInterpolator and
>> vtkExtractTimeSteps as follows, and you can set an arbitrary time
>> value to "specifiedTime" within the code:
>>
>> ----
>> #!/usr/bin/env python
>> import os
>> import vtk
>> import numpy as np
>>
>> specifiedTime = 0.18
>>
>> # Input file and variable
>> filename = os.path.abspath('mug.e')
>> nodal_var = 'convected'
>>
>> # Read Exodus Data
>> reader = vtk.vtkExodusIIReader()
>> reader.SetFileName(filename)
>> reader.UpdateInformation()
>> reader.SetTimeStep(10)
>> reader.SetAllArrayStatus(vtk.vtkExodusIIReader.NODAL, 1)
>> #reader.Update(); print reader
>>
>> info =
>> reader.GetExecutive().GetOutputInformation().GetInformationObject(0)
>> key = vtk.vtkStreamingDemandDrivenPipeline.TIME_STEPS()
>> times = np.array([info.Get(key,i) for i in range(info.Length(key))])
>> index = np.max(np.where(times <= specifiedTime))
>>
>> extractTime = vtk.vtkExtractTimeSteps()
>> extractTime.SetInputConnection(0, reader.GetOutputPort(0))
>> extractTime.SetTimeStepIndices(2, [index, index+1])
>>
>> # Time interpolation (How do I set this up?)
>> time = vtk.vtkTemporalInterpolator()
>> time.SetInputConnection(0, extractTime.GetOutputPort(0))
>> time.SetDiscreteTimeStepInterval(specifiedTime - times[index])
>> time.UpdateTimeStep(specifiedTime)
>>
>> extractTime2 = vtk.vtkExtractTimeSteps()
>> extractTime2.SetInputConnection(0, time.GetOutputPort(0))
>> extractTime2.AddTimeStepIndex(1)
>>
>> # Create Geometry
>> geometry = vtk.vtkCompositeDataGeometryFilter()
>> geometry.SetInputConnection(0, extractTime2.GetOutputPort(0))
>> geometry.Update()
>>
>> writer = vtk.vtkExodusIIWriter()
>> writer.SetInputConnection(0, geometry.GetOutputPort())
>> writer.SetFileName('time.e')
>> writer.WriteAllTimeStepsOn()
>> writer.Write()
>>
>> # Mapper
>> mapper = vtk.vtkPolyDataMapper()
>> mapper.SetInputConnection(geometry.GetOutputPort())
>> mapper.SelectColorArray(nodal_var)
>> mapper.SetScalarModeToUsePointFieldData()
>> mapper.InterpolateScalarsBeforeMappingOn()
>>
>> # Actor
>> actor = vtk.vtkActor()
>> actor.SetMapper(mapper)
>>
>> # Renderer
>> renderer = vtk.vtkRenderer()
>> renderer.AddViewProp(actor)
>>
>> # Window and Interactor
>> window = vtk.vtkRenderWindow()
>> window.AddRenderer(renderer)
>> window.SetSize(600, 600)
>>
>> interactor = vtk.vtkRenderWindowInteractor()
>> interactor.SetRenderWindow(window)
>> interactor.Initialize()
>>
>> # Show the result
>> window.Render()
>> interactor.Start()
>> ----
>>
>> Best
>>
>> 2018-05-11 2:06 GMT+09:00 Slaughter, Andrew E <andrew.slaughter at inl.gov>:
>> > I was hoping someone could help me setup interpolation for ExodusII
>> > data.
>> > The attached python scripts (I am using MacOS; python2.7; vtk7.1) reads
>> > the
>> > ExodusII file (mug.e). I would like to be able to set an arbitrary time
>> > value and have VTK perform the interpolation, would someone please help
>> > me
>> > get this working? I have two versions I am working on, one using
>> > vtkTemporalInterpolator and the other with
>> > vtkInterpolateDataSetAttributes.
>> > I am happy using either class, I just need it to work.
>> >
>> > I am aware that some of these features are built into Paraview, but my
>> > project requires me to make this work with pure VTK.
>> >
>> > Thanks,
>> >
>> > Andrew
>> >
>> > _______________________________________________
>> > Powered by www.kitware.com
>> >
>> > Visit other Kitware open-source projects at
>> >
>> > https://urldefense.proofpoint.com/v2/url?u=http-3A__www.kitware.com_opensource_opensource.html&d=DwIBaQ&c=54IZrppPQZKX9mLzcGdPfFD1hxrcB__aEkJFOKJFd00&r=h7heP8xwI1i_HikChvhFbEBurKirgfOCdwgBxB9lM8c&m=puHckBK1HPjYo6fLuDnW7Rsr9yj_a2niN1Z4KRui5BM&s=E7mpGidBh-HqQHk8zv6FqPKcN2BAfHg5FyshFS96aAM&e=
>> >
>> > Please keep messages on-topic and check the VTK FAQ at:
>> >
>> > https://urldefense.proofpoint.com/v2/url?u=http-3A__www.vtk.org_Wiki_VTK-5FFAQ&d=DwIBaQ&c=54IZrppPQZKX9mLzcGdPfFD1hxrcB__aEkJFOKJFd00&r=h7heP8xwI1i_HikChvhFbEBurKirgfOCdwgBxB9lM8c&m=puHckBK1HPjYo6fLuDnW7Rsr9yj_a2niN1Z4KRui5BM&s=__QxFPfkRQHlaX87BPaTQqC3-xiuwKePex8T9KtAD8Q&e=
>> >
>> > Search the list archives at:
>> > https://urldefense.proofpoint.com/v2/url?u=http-3A__markmail.org_search_-3Fq-3Dvtkusers&d=DwIBaQ&c=54IZrppPQZKX9mLzcGdPfFD1hxrcB__aEkJFOKJFd00&r=h7heP8xwI1i_HikChvhFbEBurKirgfOCdwgBxB9lM8c&m=puHckBK1HPjYo6fLuDnW7Rsr9yj_a2niN1Z4KRui5BM&s=CNEA4K6JVEFn6AGYM0fYeW0rg6ENnWlIx_cTPpb7zHM&e=
>> >
>> > Follow this link to subscribe/unsubscribe:
>> >
>> > https://urldefense.proofpoint.com/v2/url?u=https-3A__vtk.org_mailman_listinfo_vtkusers&d=DwIBaQ&c=54IZrppPQZKX9mLzcGdPfFD1hxrcB__aEkJFOKJFd00&r=h7heP8xwI1i_HikChvhFbEBurKirgfOCdwgBxB9lM8c&m=puHckBK1HPjYo6fLuDnW7Rsr9yj_a2niN1Z4KRui5BM&s=DJXa7N2JNLCUAVjyuFKhbN5PDsErajyMP1cjw0MJQyg&e=
>> >
>
>


More information about the vtkusers mailing list