[Paraview] questions concerning vtkFloatArray
Moreland, Kenneth
kmorel at sandia.gov
Thu Sep 11 14:41:21 EDT 2008
A couple of further comments.
The use of vtkDataSetToStructuredGridFilter is deprecated. You really
should be subclassing vtkStructuredGridAlgorithm. You override
FillInputPortInformation to report what types of input you want and override
RequestData to implement the action of the filter. The code would look
something like this:
int vtkMyFilter::FillInputPortInformation(int vtkNotUsed(port),
vtkInformation *info)
{
info->Remove(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE());
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkStructuredPoints");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkStructuredGrid");
info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(),
"vtkRectilinearGrid");
return 1;
}
int vtkMyFilter::RequestData(vtkInformation *vtkNotUsed(request),
vtkInformationVector **inputVector,
vtkInformationVector *outputVector)
{
vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
vtkStructuredGrid *output = vtkStructuredGrid::GetData(outputVector);
...
Also, you might want to consider using vtkSmartPointer when you create your
arrays and other VTK objects. That way you don't have to think about
deleting them.
-Ken
On 9/11/08 7:46 AM, "Burlen Loring" <burlen.loring at kitware.com> wrote:
> Hi Natalie,
> Some comments follow:
>
>> int dims[3],i,j,k;
>> // pointer for casting
>> vtkRectilinearGrid *rect_input = vtkRectilinearGrid::New();
>> vtkStructuredGrid *structgrid_input = vtkStructuredGrid::New();
>>
>>
>> vtkDataArray *scaldat = this ->
>> vtkAlgorithm::GetInputArrayToProcess(0,input);
>> vtkDataArray *new_scaldat;
>>
>> vtkFloatArray *scalars = vtkFloatArray::New();
>> vtkFloatArray *new_scalars = vtkFloatArray::New();
>>
>> scalars = vtkFloatArray::SafeDownCast(scaldat); //this is now the
>> input as vtkFloatArray
>
> This is going to create a memory leak. You New()'d a vtkFloatArray and
> never called Delete() before reassigning.
>
>> if((DataSetType == VTK_RECTILINEAR_GRID) && (rect_input =
>> vtkRectilinearGrid::SafeDownCast(input)))
>> {
>> rect_input -> GetDimensions(dims);
>> }
>>
>> else if((DataSetType == VTK_STRUCTURED_GRID )&& (structgrid_input =
>> vtkStructuredGrid::SafeDownCast(input)))
>> {
>> structgrid_input->GetDimensions(dims);
>> }
>>
>
> Also these are going to create leaks. You never call Delete() on the
> data sets you New()'d before reassigning.
>
>> My first problem occurs in the part, where I already got the input
>> Scalar array and casted it to vtkFloatArray. In order to perform
>> arithmetic operations with the data in it, I want to assign each point
>> of the vtkFloatArray to a "real" float array, but to do that, I need
>> to get the vtkIDList of the vtkFloatArray , so I can use the
>> GetPoint(vtkId x)-method.
>>
>> How do I get a vtkIDList of my vtkFloatArray?
> There is no GetPoint() method for a vtkFloatArray. Do you mean GetTuple?
> The Tuples are ordered such the tuple 0 is associated with point/cell 0
> for point/cell data arrays. Are you looking to get the points associated
> with the data set or the scalar values associated with the points?
>>
>> secondly,
>> i don´t really now how to set the output correctly.
>> I want to use a method like CopyStructure() so that the grid is copied
>> from input to output, but as the input is a DataSet and the output
>> StructuredGrid, does this method compute points and so for the
>> structuredGrid Format or do I have to assign them manually?
> CopyStructure, does copy the geometry. Both the source and destination
> data set have to be the same type.
>
>>
>> Thirdly, I don´t know how to assign my new vtkFloatArray as Attribute
>> Data to the output. SetScalars() does not work here, the compiler says
>> that this is not a member function of vtkStructuredGrid. But surely
>> there must be a method to assign attribute data to a structured grid?
> if you have point data then you can:
>
> output->GetPointData()->AddArray(someArray);
>
> There are similar calls for CellData and FieldData.
>
> I would advise taking a close look at the examples in the vtk user
> guide. In the filters I have seen you will implement RequestData and get
> the input and output from an instance of vtkInformation. here is an
> example using unstructured grids:
>
> int MyPv3Filter::RequestData(
> vtkInformation *req,
> vtkInformationVector **input,
> vtkInformationVector *output)
> {
> vtkDataObject *tmpDataObj;
> // Get the inputs
> tmpDataObj
> = input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
> vtkUnstructuredGrid *usgIn
> = dynamic_cast<vtkUnstructuredGrid *>(tmpDataObj);
> // Get the outputs
> tmpDataObj
> = output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
> vtkUnstructuredGrid *usgOut
> = dynamic_cast<vtkUnstructuredGrid *>(tmpDataObj);
>
> // do stuff here
>
> return 1;
> }
>
> In regard to the input arrays, they can be easily obtained from
> PointData, CellData, or FieldData via a call such as
>
> usgIn->GetPointData()->GetArray(x).
>>
>> thx,
>> Natalie
>>
>> ------------------------------------------------------------------------
>> Express yourself instantly with MSN Messenger! MSN Messenger
>> <http://clk.atdmt.com/AVE/go/onm00200471ave/direct/01/>
>> ------------------------------------------------------------------------
>>
>> _______________________________________________
>> ParaView mailing list
>> ParaView at paraview.org
>> http://www.paraview.org/mailman/listinfo/paraview
>>
>
>
> --
> Burlen Loring
> Kitware, Inc.
> R&D Engineer
> 28 Corporate Drive
> Clifton Park, NY 12065-8662
> Phone: 518-371-3971 x137
>
> _______________________________________________
> ParaView mailing list
> ParaView at paraview.org
> http://www.paraview.org/mailman/listinfo/paraview
>
More information about the ParaView
mailing list