[Paraview] [vtkusers] vtkDistributedDataFilter and ghost cells - PROGRESS

Andrew Parker andy.john.parker at googlemail.com
Wed Nov 7 11:57:14 EST 2012


Ok thanks again for the advice.  Think it's really clear I'm confused
between the "downstream" nomenclature.  Do mean just a class that
constructs or takes a D3 as input that I write from scratch outside of vtk?
Or do you mean something totally different?  Clearly my understanding of
"it was to inherit" was totally wrong.  I would prefer to stick with the
downstream filter if that's more robust.

Could you give me a quick explanation in this context?  In addition could
you point me towards an example of using a downstream filter in this
context?

Sorry again for the confusion,
Andy

On 7 November 2012 16:48, Moreland, Kenneth <kmorel at sandia.gov> wrote:

>   The instructions I gave were for a downstream filter separate from D3.
>  If you want to modify D3 itself (which sounds like a reasonable good idea
> for your purposes) then you will have change the value for
> UPDATE_NUMBER_OF_GHOST_LEVELS that the RequestData reads.  If you want to
> go with this approach, first blow away your RequestUpdateExtent method.
>  Let D3 do its regular thing.  Then implement a RequestData method like the
> following:
>
>  int MyD3::RequestData(
>    vtkInformation *request,
>   vtkInformationVector **inputVector,
>   vtkInformationVector *outputVector)
> {
>    vtkInformation *outInfo = outputVector->GetInformationObject(0);
>
>    int requestedGhostLevel = outInfo->Get(
>      vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS());
>    outInfo->Set(
>      vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS(),
>      requestedGhostLevel + 1);
>
>    return this->Superclass::RequestData(request, inputVector,
> outputVector);
> }
>
>   I don't know if the executive is going to get upset about you changing
> UPDATE_NUMBER_OF_GHOST_LEVELS in RequestData, but I think it will be all
> right.
>
>  -Ken
>
>   From: Andrew Parker <andy.john.parker at googlemail.com>
> Date: Wednesday, November 7, 2012 9:31 AM
> To: George Zagaris <george.zagaris at kitware.com>
> Cc: Kenneth Moreland <kmorel at sandia.gov>, "vtkusers at vtk.org" <
> vtkusers at vtk.org>, "paraview at paraview.org" <paraview at paraview.org>
> Subject: [EXTERNAL] Re: [vtkusers] [Paraview] vtkDistributedDataFilter
> and ghost cells - PROGRESS
>
>   So I've added this in an .cpp file and complied it:
>
>  vtkStandardNewMacro(MyD3)
>
>  MyD3::MyD3()
> : vtkDistributedDataFilter()
> {
> }
>
>  void MyD3::PrintSelf(ostream& os, vtkIndent indent)
> {
>   this->Superclass::PrintSelf(os,indent);
> }
>
>  int MyD3::RequestUpdateExtent(vtkInformation *vtkNotUsed(request),
>      vtkInformationVector **inputVector,
>      vtkInformationVector *outputVector)
> {
>   // get the info objects
>   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
>   vtkInformation *outInfo = outputVector->GetInformationObject(0);
>
>   int piece, numPieces, ghostLevels;
>
>   // We require an extra layer of ghost cells from upstream.
>
>   piece =
> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER());
>   numPieces =
> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES());
>   ghostLevels =
> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS())
> + 1; //##### New bit, modified from Ken's notes
>
>   inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER(),
> piece);
>   inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(),
> numPieces);
>
> inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS(),
> ghostLevels);
>
>   return 1;
> }
>
>  And I add this to a .h file:
>
>  class MyD3
>   : public vtkDistributedDataFilter
> {
>   vtkTypeMacro(MyD3, vtkDistributedDataFilter);
>
> public:
>   static MyD3 *New();
>   void PrintSelf(ostream& os, vtkIndent indent);
>   int RequestUpdateExtent(vtkInformation *vtkNotUsed(request),
>  vtkInformationVector **inputVector,
>  vtkInformationVector *outputVector);
> protected:
>   MyD3();
>
> private:
>   // MyD3(const MyD3&);  // Not implemented.
>   // void operator=(const MyD3&);  // Not implemented.
> };
>
>  I then invoke it in parallel 3 way, as I did before, but now using:
>
>    vtkSmartPointer<MyD3> dd = vtkSmartPointer<MyD3>::New();
>
>   // here mesh is a vtkUnstructuredGrid on each processor, obtained from
> other (non-vtk) means, but represents a valid and checked subset of the
> entire mesh.
>   dd->SetInput(mesh);
>   dd->SetController(getVtkController());
>   dd->SetBoundaryModeToSplitBoundaryCells();  // clipping
>   dd->UseMinimalMemoryOff();
>   dd->Update();
>
>   vtkSmartPointer<vtkUnstructuredGrid>
> ddUmesh(vtkUnstructuredGrid::SafeDownCast(dd->GetOutput()));
>
>   auto ids = ddUmesh->GetPointData()->GetGlobalIds();
>   auto cids = ddUmesh->GetCellData()->GetGlobalIds();
>
>  ########################
> Notes:
>
>  1) Both ids, and cids remain null, although I now understand from Ken
> that cids will always be null as they're not created.
> 2) Using totalview and stopping the code in parallel 3 way on line 670 of
> vtkDistributedDataFilter.cxx, the member this->GhostLevel is always zero
> and never executes codes branches dependent on it, however, the above code
> successfully increments the ghostLevels to one.
> 3) So basically moved forward a little bit, but still not working.  Any
> thoughts??
> 4) If I had created the ghost cells, how would I know, is my next
> question.  How can I visualise them to check?
>
>  Any thoughts really welcome.
>
>  Cheers again,
> Andy
>
>
>
> On 7 November 2012 14:50, George Zagaris <george.zagaris at kitware.com>wrote:
>
>> On Wed, Nov 7, 2012 at 9:34 AM, Andrew Parker
>> <andy.john.parker at googlemail.com> wrote:
>> > Many thanks for that.  I perform some threshold stuff early on, so
>> basically
>> >
>>
>>  Perhaps an easy solution is to append a "GlobalCellID" and
>> "GlobalNodeID" before thresholding. To the best of my knowledge,
>> threshold will extract vtkUnstructuredGrid instance from whatever goes
>> in the input and not automatically retain this type of global
>> information.
>>
>> Hope this helps.
>>
>> Best,
>> George
>>
>> > On 7 November 2012 13:00, George Zagaris <george.zagaris at kitware.com>
>> wrote:
>> >>
>> >> Hello,
>> >>
>> >> Since the data is essentially structured (rectilinear grids), I think
>> >> the best way to get the global/local mapping is to use extent
>> >> information. Specifically, you need the whole extent of the grid and
>> >> the global extent of each partition. Then from the global IJK of a
>> >> point or cell, and the whole extent you can use methods from
>> >> vtkStructuredData to compute the corresponding global linear index for
>> >> a point or cell. Do you have extent information?
>> >>
>> >> Also, I have done some work on ghost generation in VTK for structured
>> >> (and AMR) data. Unfortunately, currently it doesn't work outside the
>> >> box with rectilinear grids. This functionality can still be used to
>> >> create ghost layers and exchange the fields (node and/or cell data) on
>> >> the grid but, it doesn't communicate the coordindates on the
>> >> rectilinear grid. Alternatively, you may be able to converge to a
>> >> solution faster if you convert the rectilinear grid to a
>> >> vtkStructuredGrid and use vtkPStructuredGridConnectivity or
>> >> vtkPStructuredGridGhostDataGenerator to generate ghost layers. There
>> >> are tests that illustrate how to use these classes and I can provide
>> >> more information if you are interested.
>> >>
>> >> Best,
>> >> George
>> >>
>> >> On Wed, Nov 7, 2012 at 7:22 AM, Moreland, Kenneth <kmorel at sandia.gov>
>> >> wrote:
>> >> > Hmmm. I thought global point IDs were always created. Maybe it only
>> >> > happens
>> >> > when ghost cells are created. Try getting that working first.
>> >> >
>> >> > As far as global cell IDs are concerned, you are out of luck.
>> >> > vtkDistributedDataFilter does not create them. But they are not all
>> that
>> >> > useful anyway since, unlike points, once you remove the out cells
>> there
>> >> > is
>> >> > no overlap. Usually you can get all the information you want from the
>> >> > ghost
>> >> > level.
>> >> >
>> >> > If you really need global cell IDs, I'm sure you would have no
>> trouble
>> >> > making them yourself. I would do it before running
>> >> > vtkDistributedDataFilter.
>> >> > Assuming when you first read in your data you have no duplicate
>> cells,
>> >> > just
>> >> > label those on process 0 as 0..N_0, those on process 1 as N_0+1..N_1,
>> >> > and so
>> >> > on.
>> >> >
>> >> > -Ken
>> >> >
>> >> > Sent from my iPad so blame autocorrect.
>> >> >
>> >> > On Nov 7, 2012, at 2:38 AM, "Andrew Parker"
>> >> > <andy.john.parker at googlemail.com> wrote:
>> >> >
>> >> > So that's working now in terms of the cast, forgot GetOutput() inside
>> >> > the
>> >> > cast operator!  The returned vtkUGrid is now fully populated and
>> >> > contains
>> >> > sensible information.
>> >> >
>> >> > However, both getPointData()->GetGlobalIds()  and
>> >> > getCellData()->GetGlobalIds() return null pointers.  Any thoughts?
>> >> >
>> >> > Also, should I be using CellData since I want the cell global to
>> local
>> >> > mapping for cells not the nodes, at the moment at least?
>> >> >
>> >> > On 6 November 2012 19:05, Moreland, Kenneth <kmorel at sandia.gov>
>> wrote:
>> >> >>
>> >> >> Perhaps it is outputting a composite data set of some type?  Try
>> >> >> running
>> >> >> GetClassName() to see what type of data object it really is.
>> >> >>
>> >> >> -Ken
>> >> >>
>> >> >> From: Andrew Parker <andy.john.parker at googlemail.com>
>> >> >> Date: Tuesday, November 6, 2012 9:50 AM
>> >> >>
>> >> >> To: Kenneth Moreland <kmorel at sandia.gov>
>> >> >> Cc: "vtkusers at vtk.org" <vtkusers at vtk.org>, "paraview at paraview.org"
>> >> >> <paraview at paraview.org>
>> >> >> Subject: [EXTERNAL] Re: [Paraview] vtkDistributedDataFilter and
>> ghost
>> >> >> cells
>> >> >>
>> >> >> Thanks on both accounts.  Any thoughts why the downcast called after
>> >> >> dd->Update() on distributedDataFilter is a null pointer?  As in, dd
>> is
>> >> >> working perfectly properly, but I don't seem to be able to extract a
>> >> >> valid
>> >> >> unstructuredgrid.  For a follow up question I assume
>> >> >> getPointData()->GetGlobalIds() gives the local to global for the
>> mesh
>> >> >> nodes,
>> >> >> and I should use getCellData()->GetGlobalIds() to get the local to
>> >> >> global
>> >> >> for the cells?  Once I get a valid pointer that is....
>> >> >>
>> >> >> Cheers again,
>> >> >> Andy
>> >> >>
>> >> >> On 6 November 2012 16:40, Moreland, Kenneth <kmorel at sandia.gov>
>> wrote:
>> >> >>>
>> >> >>> You should be able to do a vtkUnstructuredGrid::SafeDownCast() to
>> the
>> >> >>> data to get the unstructured mesh (and access to the point data).
>> >> >>>
>> >> >>> -Ken
>> >> >>>
>> >> >>> From: Andrew Parker <andy.john.parker at googlemail.com>
>> >> >>> Date: Tuesday, November 6, 2012 9:32 AM
>> >> >>> To: Kenneth Moreland <kmorel at sandia.gov>
>> >> >>> Cc: "vtkusers at vtk.org" <vtkusers at vtk.org>, "paraview at paraview.org"
>> >> >>> <paraview at paraview.org>
>> >> >>> Subject: [EXTERNAL] Re: [Paraview] vtkDistributedDataFilter and
>> ghost
>> >> >>> cells
>> >> >>>
>> >> >>> Another question which I'd forgotten, how do I get to a
>> >> >>> vtkUnstructuredGrid per processor from the
>> vtkDistributedDataFilter.
>> >> >>>
>> >> >>> For instance, dd->GetOutput()->GetPointData()->GetGlobalIds()
>> doesn't
>> >> >>> work as it's a vtkDataObject
>> >> >>>
>> >> >>> Stupid question I'm sure but the doxy notes say this type returns
>> an
>> >> >>> unstructured mesh, but I can't seem to get it out?
>> >> >>>
>> >> >>> Also, why exactly do I need the vtkPieceScalars and
>> >> >>> vtkDataSetSurfaceFilter again? If the above can be made to work and
>> >> >>> return
>> >> >>> the mapping, what are they adding in terms of information?
>> >> >>>
>> >> >>> Thanks again,
>> >> >>> Andy
>> >> >>>
>> >> >>> On 6 November 2012 16:00, Moreland, Kenneth <kmorel at sandia.gov>
>> wrote:
>> >> >>>>
>> >> >>>> I believe vtkDistributedDataFilter will always return with global
>> >> >>>> point
>> >> >>>> ids (a mapping from local point ids to global point ids),
>> although it
>> >> >>>> might
>> >> >>>> pass them if they already exist.  So
>> >> >>>> dd->GetOutput()->GetPointData()->GetGlobalIds() should return the
>> >> >>>> array that
>> >> >>>> gives this mapping.
>> >> >>>>
>> >> >>>> Ghost cells are only created on demand, and this is usually done
>> by
>> >> >>>> pipeline convention.  If you have a filter that needs a layer of
>> >> >>>> ghost
>> >> >>>> cells, it should override the RequestUpdateExtent method to
>> increment
>> >> >>>> the
>> >> >>>> value of UPDATE_NUMBER_OF_GHOST_LEVELS from the output
>> information to
>> >> >>>> the
>> >> >>>> input information.  This method would look something like this.
>> >> >>>>
>> >> >>>> int vtkDistributedDataFilter::RequestUpdateExtent(
>> >> >>>>   vtkInformation *vtkNotUsed(request),
>> >> >>>>   vtkInformationVector **inputVector,
>> >> >>>>   vtkInformationVector *outputVector)
>> >> >>>> {
>> >> >>>>   // get the info objects
>> >> >>>>   vtkInformation *inInfo =
>> inputVector[0]->GetInformationObject(0);
>> >> >>>>   vtkInformation *outInfo = outputVector->GetInformationObject(0);
>> >> >>>>
>> >> >>>>   int piece, numPieces, ghostLevels;
>> >> >>>>
>> >> >>>>   // We require an extra layer of ghost cells from upstream.
>> >> >>>>
>> >> >>>>   piece =
>> >> >>>>
>> >> >>>>
>> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER());
>> >> >>>>   numPieces =
>> >> >>>>
>> >> >>>>
>> >> >>>>
>> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES());
>> >> >>>>   ghostLevels =
>> >> >>>>
>> >> >>>>
>> >> >>>>
>> outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS());
>> >> >>>>
>> >> >>>>
>> >> >>>>
>> inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER(),
>> >> >>>> piece);
>> >> >>>>
>> >> >>>>
>> >> >>>>
>> inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(),
>> >> >>>>               numPieces);
>> >> >>>>
>> >> >>>>
>> >> >>>>
>> inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_GHOST_LEVELS(),
>> >> >>>>               ghostLevels);
>> >> >>>>
>> >> >>>>   return 1;
>> >> >>>> }
>> >> >>>>
>> >> >>>>
>> >> >>>> The operation of the RequestData method should also strip off this
>> >> >>>> layer
>> >> >>>> of ghost cells.  It might be possible to request a layer of ghost
>> >> >>>> cells by
>> >> >>>> setting UPDATE_NUMBER_OF_GHOST_LEVELS at the bottom of the
>> pipeline,
>> >> >>>> but I'm
>> >> >>>> not totally sure how to make that work.  It's probably easier (or
>> at
>> >> >>>> least
>> >> >>>> cleaner) to do it from within a filter.
>> >> >>>>
>> >> >>>> -Ken
>> >> >>>>
>> >> >>>> From: Andrew Parker <andy.john.parker at googlemail.com>
>> >> >>>> Date: Tuesday, November 6, 2012 8:25 AM
>> >> >>>> To: "vtkusers at vtk.org" <vtkusers at vtk.org>, "paraview at paraview.org
>> "
>> >> >>>> <paraview at paraview.org>
>> >> >>>> Subject: [EXTERNAL] [Paraview] vtkDistributedDataFilter and ghost
>> >> >>>> cells
>> >> >>>>
>> >> >>>> Hi,
>> >> >>>>
>> >> >>>> Hope you can help.  I have some code running in parallel, that by
>> >> >>>> other
>> >> >>>> means I have constructed nprocs worth of vtkRectilinearGrids, one
>> per
>> >> >>>> process.  Each of which is a valid nprocs-worth of the whole
>> serial
>> >> >>>> mesh,
>> >> >>>> I've check this and I am happy with that i.e. it's partitioned
>> >> >>>> properly and
>> >> >>>> nothing is missing.  I need the following information to process
>> my
>> >> >>>> data in
>> >> >>>> parallel:
>> >> >>>>
>> >> >>>> 1) I would like the local -> global cell mapping between the local
>> >> >>>> rgrid
>> >> >>>> and the corresponding global single mesh.
>> >> >>>> 2) I would like to know which cells are on processor boundaries
>> for
>> >> >>>> parallel exchange purposes.
>> >> >>>> 3) I would like all the double arrays per processor to be
>> "expanded"
>> >> >>>> by
>> >> >>>> the amount of (1 level of) ghost cells such that I can properly do
>> >> >>>> the
>> >> >>>> computations I want with the ability to exchange only those
>> >> >>>> additional cells
>> >> >>>> given the local to global mapping.
>> >> >>>>
>> >> >>>> I have tried from the examples to use the following code, which I
>> >> >>>> call
>> >> >>>> on every process, each of which has it's own local rgrid as I
>> said.
>> >> >>>> I do
>> >> >>>> the following:
>> >> >>>>
>> >> >>>>  vtkSmartPointer<vtkDistributedDataFilter> dd =
>> >> >>>> vtkSmartPointer<vtkDistributedDataFilter>::New();
>> >> >>>>  dd->SetInput(rgrid);
>> >> >>>>
>> >> >>>>  dd->SetController(getVtkController());
>> >> >>>>  dd->SetBoundaryModeToSplitBoundaryCells();
>> >> >>>>  //dd->SetBoundaryModeToAssignToOneRegion();
>> >> >>>>  //dd->SetBoundaryModeToAssignToAllIntersectingRegions();
>> >> >>>>  dd->UseMinimalMemoryOff();
>> >> >>>>  dd->Update();
>> >> >>>>  vtkPieceScalars *ps = vtkPieceScalars::New();
>> >> >>>>  ps->SetInputConnection(dd->GetOutputPort());
>> >> >>>>  ps->SetScalarModeToCellData();
>> >> >>>>  vtkDataSetSurfaceFilter *dss = vtkDataSetSurfaceFilter::New();
>> >> >>>>  dss->SetInputConnection(ps->GetOutputPort());
>> >> >>>>
>> >> >>>> The dd object works fine and writing its contents out on each
>> >> >>>> processor
>> >> >>>> gives nprocs worth of meshes, each of which look slightly
>> different
>> >> >>>> to the
>> >> >>>> way I've partitioned them up, but sum to the same serial mesh so
>> I am
>> >> >>>> happy
>> >> >>>> with that working correctly. But I can't for the life of me figure
>> >> >>>> out how
>> >> >>>> to obtain local to global cell mappings, allocate ghost cells, or
>> >> >>>> work out
>> >> >>>> how to exchange data given the above partition info and comms....
>> >> >>>>
>> >> >>>> Note I have not provided any additional information to "dd"
>> regarding
>> >> >>>> global cells as per the doxy notes so I assume it went away and
>> >> >>>> computed it.
>> >> >>>> I can't figure out how to extract it however.  I also have no idea
>> >> >>>> how to
>> >> >>>> modify each local processor rgrid to include the ghost cells for
>> that
>> >> >>>> processor.  Finally given that info, I could exchange between
>> >> >>>> processors to
>> >> >>>> write to each local processors ghost cells the corresponding
>> "real"
>> >> >>>> cell
>> >> >>>> data from the neighbouring meshes and continue the code.
>> >> >>>>
>> >> >>>> Any help really appreciated!
>> >> >>>>
>> >> >>>> Cheers,
>> >> >>>> Andy
>> >> >>>>
>> >> >>>
>> >> >>>
>> >> >>>
>> >> >>> --
>> >> >>>
>> >> >>> __________________________________
>> >> >>>
>> >> >>>    Dr Andrew Parker
>> >> >>>
>> >> >>>    Em at il:  andrew.parker at cantab.net
>> >> >>>
>> >> >>>
>> >> >>
>> >> >>
>> >> >>
>> >> >> --
>> >> >>
>> >> >> __________________________________
>> >> >>
>> >> >>    Dr Andrew Parker
>> >> >>
>> >> >>    Em at il:  andrew.parker at cantab.net
>> >> >>
>> >> >>
>> >> >
>> >> >
>> >> >
>> >> > --
>> >> >
>> >> > __________________________________
>> >> >
>> >> >    Dr Andrew Parker
>> >> >
>> >> >    Em at il:  andrew.parker at cantab.net
>> >> >
>> >> >
>> >> >
>> >> > _______________________________________________
>> >> > Powered by www.kitware.com
>> >> >
>> >> > Visit other Kitware open-source projects at
>> >> > http://www.kitware.com/opensource/opensource.html
>> >> >
>> >> > Please keep messages on-topic and check the VTK FAQ at:
>> >> > http://www.vtk.org/Wiki/VTK_FAQ
>> >> >
>> >> > Follow this link to subscribe/unsubscribe:
>> >> > http://www.vtk.org/mailman/listinfo/vtkusers
>> >> >
>> >
>> >
>> >
>> >
>> > --
>> >
>> > __________________________________
>> >
>> >    Dr Andrew Parker
>> >
>> >    Em at il:  andrew.parker at cantab.net
>> >
>> >
>>
>
>
>
>  --
>
> __________________________________
>
>    Dr Andrew Parker
>
>    Em at il:  andrew.parker at cantab.net
>
>


-- 

__________________________________

   Dr Andrew Parker

   Em at il:  andrew.parker at cantab.net
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://www.paraview.org/pipermail/paraview/attachments/20121107/0ce5c008/attachment-0001.htm>


More information about the ParaView mailing list