Multi-Resolution Rendering with Overlapping AMR

From ParaQ Wiki
Jump to navigationJump to search

Motivation

With dataset sizes growing increasingly large, exploratory interactive visualization keeps on getting harder and harder. Multi-resoltuion visualization makes it possible to deal with large datasets while keeping the system requirements and response times low.

This document is for advanced developers and assumes familiarity with VTK data and execution model and ParaView architecture.

Design Overview

ParaView/VTK provides an infrastructure for multi-resolution analysis based on composite datasets (multiblock datasets or overlapping AMR datasets). The principle is as follows:

  1. Sources/Filters produce meta-data in their RequestInformation() pass. This meta-data can be used by any sink to decide what are the blocks in the dataset, how are they positioned, what is the cost for processing that block of dataset. For AMR datasets, this meta-data is properly defined by vtkAMRInformation. For other composite datasets, users are free to come up with their own convention. So long as the sink and source are aware of the convention, we are good.
  2. A Sink can request a particular block (or blocks) from the source by providing appropriate keys in the RequestUpdateExtent() pass. The source then delivers those blocks.

This general streaming principal can be extended to views and representations too, for rendering. The representations, in this case, act as the sinks that process meta-data and request blocks from the sources. The render view in ParaView (vtkPVRenderView), provides mechanisms for representations to leverage to render the dataset in a streaming fashion. In additional to the usual view passes, when a representation indicates that it supports streaming (and streaming is enabled), the view adds two additional passes vtkPVRenderView::REQUEST_STREAMING_UPDATE() in which a representation can request the "next" block, possibly using the current camera positions, and vtkPVRenderView::REQUEST_PROCESS_STREAMED_PIECE() in which a representation can "handle" the streaming block to put it into the rendering pipeline. The view takes over control between REQUEST_STREAMING_UPDATE() and REQUEST_PROCESS_STREAMED_PIECE() to ensure that the data is delivered to the rendering processes correctly.

Implementation

While the framework is flexible enough to support any composite dataset, our current implementation only supports AMR datasets (vtkOverlappingAMR). This section provides an overview of the various components involved in multi-resolution rendering of AMR datasets.

Data Producer

The starting point is the data source. The source needs to provide meta-data about the data in RequestInformation() pass and then read blocks of data as requested in the RequestData() pass. Examples are vtkUniformGridAMRReader, and vtkAMRFlashReader. These readers read the "meta-data" in the files (without reading the heavy data) and produce an empty vtkOverlappingAMR in their RequestInformation pass with structure matching the expected blocks in the output.

//----------------------------------------------------------------------------
int vtkXMLUniformGridAMRReader::RequestInformation(vtkInformation *request,
  vtkInformationVector **inputVector, vtkInformationVector *outputVector)
{
  if (!this->Superclass::RequestInformation(request, inputVector, outputVector))
    {
    return 0;
    }

  // this->Metadata is of type vtkOverlappingAMR and is filled up before the control reaches
  // this point.
  if (this->Metadata)
    {
    // Pass the meta-data down the pipeline.
    outputVector->GetInformationObject(0)->Set(
      vtkCompositeDataPipeline::COMPOSITE_DATA_META_DATA(),
      this->Metadata);
    }
  else
    {
    outputVector->GetInformationObject(0)->Remove(
      vtkCompositeDataPipeline::COMPOSITE_DATA_META_DATA());
    }
  return 1;
}

Data Consumer / Representation

Now we need to provide a representation that can use this meta-data to request blocks in a streaming fashion. This representation can have any arbitrary code to manage multiple blocks and thus can employ arbitrary conventions about how the data is structured, how multiple resolutions are merged, etc. Examples are vtkAMROutlineRepresentation, and vtkAMRStreamingVolumeRepresentation. vtkAMROutlineRepresentation is a simple representation that simply renders blocks for datasets delivered. It is a good illustration for representation that support streaming since the rendering code itself is simple.

Let's look at vtkAMROutlineRepresentation in more detail.

First, this representation needs to check that the input pipeline is streaming capable. This is done by overriding RequestInformation() and checking for presence of the meta-data. vtkPVView::GetEnableStreaming() maybe removed in future. It is currently provided to allow enabling/disabling streaming support from command line. vtkPVView::GetEnableStreaming() returns true when the process is started with "--enable-streaming" command line argument.

//----------------------------------------------------------------------------
int vtkAMROutlineRepresentation::RequestInformation(vtkInformation *rqst,
    vtkInformationVector **inputVector,
    vtkInformationVector *outputVector)
{
  // Determine if the input is streaming capable. A pipeline is streaming
  // capable if it provides us with COMPOSITE_DATA_META_DATA() in the
  // RequestInformation() pass. It implies that we can request arbitrary blocks
  // from the input pipeline which implies stream-ability.

  this->StreamingCapablePipeline = false;
  if (inputVector[0]->GetNumberOfInformationObjects() == 1)
    {
    vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
    if (inInfo->Has(vtkCompositeDataPipeline::COMPOSITE_DATA_META_DATA()) &&
      vtkPVView::GetEnableStreaming())
      {
      this->StreamingCapablePipeline = true;
      }
    }

  vtkStreamingStatusMacro(
    << this << ": streaming capable input pipeline? "
    << (this->StreamingCapablePipeline? "yes" : "no"));
  return this->Superclass::RequestInformation(rqst, inputVector, outputVector);
}

Next, we need to tell the view that this representation supports streaming. This is done in ProcessViewRequest() as follows:

//----------------------------------------------------------------------------
int vtkAMROutlineRepresentation::ProcessViewRequest(
  vtkInformationRequestKey* request_type, vtkInformation* inInfo, vtkInformation* outInfo)
{
  // always forward to superclass first. Superclass returns 0 if the
  // representation is not visible (among other things). In which case there's
  // nothing to do.
  if (!this->Superclass::ProcessViewRequest(request_type, inInfo, outInfo))
    {
    return 0;
    }

  if (request_type == vtkPVView::REQUEST_UPDATE())
    {
    // Standard representation stuff, first.
    // 1. Provide the data being rendered.
    vtkPVRenderView::SetPiece(inInfo, this, this->ProcessedData);

    // 2. Provide the bounds.
    double bounds[6];
    this->DataBounds.GetBounds(bounds);
    vtkPVRenderView::SetGeometryBounds(inInfo, bounds);

    // The only thing extra we need to do here is that we need to let the view
    // know that this representation is streaming capable (or not).
    vtkPVRenderView::SetStreamable(inInfo, this, this->GetStreamingCapablePipeline());
    }
  
  ...
}

In RequestData(), the representation is expected to convert the input data to the representable form, outlines in our case. Since RequestData() is called in two cases when the input data changes i.e. pass 0, and additional streaming passes, this code typically has to be smart enough to realize if we are in vtkPVRenderView::REQUEST_STREAMING_UPDATE() pass. When not in REQUEST_STREAMING_UPDATE(), vtkAMROutlineRepresentation initializes internal data-structure that computes priorities from blocks using their bounds, level, and current camera parameters using the uility class vtkAMRStreamingPriorityQueue.

int vtkAMROutlineRepresentation::RequestData(vtkInformation *rqst,
  vtkInformationVector **inputVector, vtkInformationVector *outputVector)
{
  if (inputVector[0]->GetNumberOfInformationObjects() == 1)
    {
    vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
    if (inInfo->Has(vtkCompositeDataPipeline::COMPOSITE_DATA_META_DATA()) &&
      this->GetStreamingCapablePipeline() &&
      !this->GetInStreamingUpdate())
      {
      // Since the representation reexecuted, it means that the input changed
      // and we should initialize our streaming.
      vtkOverlappingAMR* amr = vtkOverlappingAMR::SafeDownCast(
        inInfo->Get(vtkCompositeDataPipeline::COMPOSITE_DATA_META_DATA()));
      this->PriorityQueue->Initialize(amr->GetAMRInfo());
      }
    }

  this->ProcessedPiece = 0;
  if (inputVector[0]->GetNumberOfInformationObjects() == 1)
    {
    // Do the streaming independent "transformation" of the data here, in our
    // case, generate the outline from the input data.
    // To keep things simple here, we don't bother about the "flip-book" caching
    // support used for animation playback.

    vtkNew<vtkPVGeometryFilter> geomFilter;
    geomFilter->SetUseOutline(1);
    geomFilter->SetHideInternalAMRFaces(false);
    // since we want to show the real data boundaries, we force the filter to
    // not use the meta-data.
    geomFilter->SetUseNonOverlappingAMRMetaDataForOutlines(false);

    vtkOverlappingAMR* input = vtkOverlappingAMR::GetData(inputVector[0], 0);
    geomFilter->SetInputData(input);
    geomFilter->Update();
    if (!this->GetInStreamingUpdate())
      {
      this->ProcessedData = geomFilter->GetOutputDataObject(0);

      double bounds[6];
      input->GetBounds(bounds);
      this->DataBounds.SetBounds(bounds);
      }
    else
      {
      this->ProcessedPiece = geomFilter->GetOutputDataObject(0);
      }
    }
  else
    {
    // create an empty dataset. This is needed so that view knows what dataset
    // to expect from the other processes on this node.
    this->ProcessedData = vtkSmartPointer<vtkMultiBlockDataSet>::New();
    this->DataBounds.Reset();
    }

  if (!this->GetInStreamingUpdate())
    {
    this->RenderedData = 0;

    // provide the mapper with an empty input. This is needed only because
    // mappers die when input is NULL, currently.
    vtkNew<vtkMultiBlockDataSet> tmp;
    this->Mapper->SetInputDataObject(tmp.GetPointer());
    }

  return this->Superclass::RequestData(rqst, inputVector, outputVector);
}

Further in ProcessViewRequest(), we need to handle the streaming passes REQUEST_STREAMING_UPDATE(), where we make a request for a block based on current view paramaters, and REQUEST_PROCESS_STREAMED_PIECE(), where we merge the block being delivered into the data we already have for rendering. This is done as follows.

//----------------------------------------------------------------------------
int vtkAMROutlineRepresentation::ProcessViewRequest(
  vtkInformationRequestKey* request_type, vtkInformation* inInfo, vtkInformation* outInfo)
{

  ....

  else if (request_type == vtkPVRenderView::REQUEST_STREAMING_UPDATE())
    {
    if (this->GetStreamingCapablePipeline())
      {
      // This is a streaming update request, request next piece.
      double view_planes[24];
      inInfo->Get(vtkPVRenderView::VIEW_PLANES(), view_planes);
      if (this->StreamingUpdate(view_planes))
        {
        // since we indeed "had" a next piece to produce, give it to the view
        // so it can deliver it to the rendering nodes.
        vtkPVRenderView::SetNextStreamedPiece(
          inInfo, this, this->ProcessedPiece);
        }
      }
    }
  else if (request_type == vtkPVRenderView::REQUEST_PROCESS_STREAMED_PIECE())
    {
    vtkDataObject* piece = vtkPVRenderView::GetCurrentStreamedPiece(inInfo, this);
    if (piece)
      {
      assert (this->RenderedData != NULL);
      vtkStreamingStatusMacro( << this << ": received new piece.");

      // merge with what we are already rendering.
      vtkNew<vtkAppendCompositeDataLeaves> appender;
      appender->AddInputDataObject(piece);
      appender->AddInputDataObject(this->RenderedData);
      appender->Update();

      this->RenderedData = appender->GetOutputDataObject(0);
      this->Mapper->SetInputDataObject(this->RenderedData);
      }
    }

  return 1;
}

vtkAMROutlineRepresentation::StreamingUpdate(), used above, is fairly simple piece of code. It updates the priority queue and triggers an update as follows.

bool vtkAMROutlineRepresentation::StreamingUpdate(const double view_planes[24])
{
  assert(this->InStreamingUpdate == false);
  if (!this->PriorityQueue->IsEmpty())
    {
    this->InStreamingUpdate = true;
    vtkStreamingStatusMacro(<< this << ": doing streaming-update.")

    // update the priority queue, if needed.
    this->PriorityQueue->Update(view_planes);

    // This ensure that the representation re-executes.
    this->MarkModified();

    // Execute the pipeline.
    this->Update();

    this->InStreamingUpdate = false;
    return true;
    }

  return false;
}

In RequestUpdateExtent(), we then request appropriate blocks if the update was called while streaming was in progress.

//----------------------------------------------------------------------------
int vtkAMROutlineRepresentation::RequestUpdateExtent(
  vtkInformation* request,
  vtkInformationVector** inputVector,
  vtkInformationVector* outputVector)
{
  if (!this->Superclass::RequestUpdateExtent(request, inputVector,
      outputVector))
    {
    return 0;
    }

  for (int cc=0; cc < this->GetNumberOfInputPorts(); cc++)
    {
    for (int kk=0; kk < inputVector[cc]->GetNumberOfInformationObjects(); kk++)
      {
      vtkInformation* info = inputVector[cc]->GetInformationObject(kk);
      if (this->InStreamingUpdate)
        {
        assert(this->PriorityQueue->IsEmpty() == false);
        int cid = static_cast<int>(this->PriorityQueue->Pop());
        vtkStreamingStatusMacro(<< this << ": requesting blocks: " << cid);
        // Request the next "group of blocks" to stream.
        info->Set(vtkCompositeDataPipeline::LOAD_REQUESTED_BLOCKS(), 1);
        info->Set(vtkCompositeDataPipeline::UPDATE_COMPOSITE_INDICES(), &cid, 1);
        }
      else
        {
        // let the source deliver whatever is the default. What the reader does
        // when the downstream doesn't request any particular blocks in poorly
        // defined right now. I am assuming the reader will only read the root
        // block or down to some user-specified level.
        info->Remove(vtkCompositeDataPipeline::LOAD_REQUESTED_BLOCKS());
        info->Remove(vtkCompositeDataPipeline::UPDATE_COMPOSITE_INDICES());
        }
      }
    }
  return 1;
}