[Paraview] Calculator speed
pat marion
pat.marion at kitware.com
Thu Sep 11 16:48:42 EDT 2008
Or in the meantime...
>>> from paraview import numpy_support
>>> help(numpy_support)
Help on module paraview.numpy_support in paraview:
NAME
paraview.numpy_support
FILE
/source/paraview/streaming-release/build/Utilities/VTKPythonWrapping/paraview/numpy_support.py
DESCRIPTION
This module adds support to easily import and export NumPy
(http://numpy.scipy.org) arrays into/out of VTK arrays. The code is
loosely based on TVTK (https://svn.enthought.com/enthought/wiki/TVTK).
This code depends on an addition to the VTK data arrays made by Berk
Geveci to make it support Python's buffer protocol (on Feb. 15, 2008).
The main functionality of this module is provided by the two functions:
numpy_to_vtk,
vtk_to_numpy.
Caveats:
--------
- Bit arrays in general do not have a numpy equivalent and are not
supported. Char arrays are also not easy to handle and might not
work as you expect. Patches welcome.
- You need to make sure you hold a reference to a Numpy array you want
to import into VTK. If not you'll get a segfault (in the best case).
The same holds in reverse when you convert a VTK array to a numpy
array -- don't delete the VTK array.
Created by Prabhu Ramachandran in Feb. 2008.
FUNCTIONS
create_vtk_array(vtk_arr_type)
Internal function used to create a VTK data array from another
VTK array given the VTK array type.
get_numpy_array_type(vtk_array_type)
Returns a numpy array typecode given a VTK array type.
get_vtk_array_type(numpy_array_type)
Returns a VTK typecode given a numpy array.
get_vtk_to_numpy_typemap()
Returns the VTK array type to numpy array type mapping.
numpy_to_vtk(num_array, deep=0)
Converts a contiguous real numpy Array to a VTK array object.
This function only works for real arrays that are contiguous.
Complex arrays are NOT handled. It also works for multi-component
arrays. However, only 1, and 2 dimensional arrays are supported.
This function is very efficient, so large arrays should not be a
problem.
If the second argument is set to 1, the array is deep-copied from
from numpy. This is not as efficient as the default behavior
(shallow copy) and uses more memory but detaches the two arrays
such that the numpy array can be released.
WARNING: You must maintain a reference to the passed numpy array, if
the numpy data is gc'd and VTK will point to garbage which will in
the best case give you a segfault.
Parameters
----------
- num_array : a contiguous 1D or 2D, real numpy array.
vtk_to_numpy(vtk_array)
Converts a VTK data array to a numpy array.
Given a subclass of vtkDataArray, this function returns an
appropriate numpy array containing the same data -- it actually
points to the same data.
WARNING: This does not work for bit arrays.
Parameters
----------
- vtk_array : `vtkDataArray`
The VTK data array to be converted.
DATA
VTK_ID_TYPE_SIZE = 8
VTK_LONG_TYPE_SIZE = 8
On Thu, Sep 11, 2008 at 2:09 PM, Berk Geveci <berk.geveci at kitware.com> wrote:
> It will be once we update the wiki pages :-)
>
> On Thu, Sep 11, 2008 at 1:18 PM, Weirs, V Gregory <vgweirs at sandia.gov> wrote:
>>
>> Dig it.
>>
>> Is numpy_support documented anywhere?
>>
>> Thanks,
>> Greg
>>
>>
>> On 9/11/08 10:31 AM, "Berk Geveci" <berk.geveci at kitware.com> wrote:
>>
>> (assuming you are using cvs paraview)
>> In the programmable filter, you can do something like:
>>
>> from paraview import numpy_support
>>
>> input = self.GetInputDataObject(0, 0)
>> ar = numpy_support.vtk_to_numpy(input.GetPointData().GetArray("pressure"))
>> ar2 = ar*2
>> newArray = nump7_support.numpy_to_vtk(ar2, True)
>> newArray.SetName("pressure times 2")
>> output = self.GetOutputDataObject(0)
>> output.GetPointData().AddArray(newArray)
>>
>> In the future, this should look like:
>> output.point_data['pressure times 2'] = input.point_data['pressure'] * 2
>>
>> -berk
>>
>>
>> On Tue, Sep 9, 2008 at 2:40 PM, Tobias Brandvik <tbrandvik at gmail.com> wrote:
>>> Thanks for the replies. I've got some experience with Python/NumPy so
>>> any pointers on that approach would be appreciated.
>>>
>>> Tobias
>>>
>>>
>>> --
>>> Tobias Brandvik
>>> PhD Student
>>> Whittle Laboratory
>>> 1 JJ Thomson Avenue
>>> Cambridge CB3 0DY, UK
>>>
>>> On Tue, Sep 9, 2008 at 2:51 PM, Berk Geveci <berk.geveci at kitware.com>
>>> wrote:
>>>> The array calculator does not use Python but it is slow nevertheless.
>>>> Actually, this could be done much faster using Python/NumPy. This is
>>>> currently possible but hard to do (complicated API). We are working on
>>>> an easier-to-use interface. If you are adventurous, I can give you
>>>> some pointers. Otherwise, C++ custom filter is the way to go.
>>>>
>>>> -berk
>>>>
>>>> On Tue, Sep 9, 2008 at 8:58 AM, Tobias Brandvik <tbrandvik at gmail.com>
>>>> wrote:
>>>>> Hi,
>>>>>
>>>>> I've written a CFD solver that uses Paraview in parallel mode with
>>>>> hdf5/xdmf. Currently, the application writes out the primary variables
>>>>> and lets the user calculate any desired secondary variables using the
>>>>> calculator. However, this procedure is very slow -- it seems to be an
>>>>> order of magnitude faster to calculate the data in the solver itself
>>>>> and have Paraview read it from disk. I suspect this is because the
>>>>> calculator expressions get translated to python code, but I might also
>>>>> be doing something wrong.
>>>>>
>>>>> Am I correct in thinking that the easiest way to fix this would be to
>>>>> write a custom filter in C++?
>>>>>
>>>>> Cheers,
>>>>> Tobias
>>>>>
>>>>> --
>>>>> Tobias Brandvik
>>>>> PhD Student
>>>>> Whittle Laboratory
>>>>> 1 JJ Thomson Avenue
>>>>> Cambridge CB3 0DY, UK
>>>>> _______________________________________________
>>>>> ParaView mailing list
>>>>> ParaView at paraview.org
>>>>> http://www.paraview.org/mailman/listinfo/paraview
>>>>>
>>>>
>>>
>> _______________________________________________
>> ParaView mailing list
>> ParaView at paraview.org
>> http://www.paraview.org/mailman/listinfo/paraview
>>
>>
>>
> _______________________________________________
> ParaView mailing list
> ParaView at paraview.org
> http://www.paraview.org/mailman/listinfo/paraview
>
More information about the ParaView
mailing list