VTK
vtkStructuredImplicitConnectivity.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkStructuredImplicitConnectivity.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
42#ifndef vtkStructuredImplicitConnectivity_h
43#define vtkStructuredImplicitConnectivity_h
44
45#include "vtkFiltersParallelMPIModule.h" // For export macro
46#include "vtkObject.h"
47
48// Forward declarations
49class vtkDataArray;
50class vtkImageData;
53class vtkPointData;
54class vtkPoints;
57
58namespace vtk
59{
60namespace detail
61{
62
63class CommunicationManager;
64struct DomainMetaData;
65struct StructuredGrid;
66
67} // END namespace detail
68} // END namespace vtk
69
70class VTKFILTERSPARALLELMPI_EXPORT vtkStructuredImplicitConnectivity :
71 public vtkObject
72{
73public:
75 void PrintSelf(ostream& os, vtkIndent indent);
77
84 void SetWholeExtent(int wholeExt[6]);
85
86 // \brief Registers the structured grid dataset belonging to this process.
87 // \param gridID the ID of the grid in this rank.
88 // \param extent the [imin,imax,jmin,jmax,kmin,kmax] of the grid.
89 // \param gridPnts pointer to the points of the grid (NULL for uniform grid).
90 // \param pointData pointer to the node-centered fields of the grid.
91 // \pre gridID >= 0. The code uses values of gridID < -1 as flag internally.
92 // \pre vtkStructuredExtent::Smaller(extent,wholeExtent) == true.
93 // \note A rank with no or an empty grid, should not call this method.
95 const int gridID,
96 int extent[6],
97 vtkPoints* gridPnts,
98 vtkPointData* pointData
99 );
100
101 // \brief Registers the rectilinear grid dataset belonging to this process.
102 // \param gridID the ID of the in this rank.
103 // \param extent the [imin,imax,jmin,jmax,kmin,kmax] of the grid.
104 // \param xcoords the x-coordinates array of the rectilinear grid.
105 // \param ycoords the y-coordinates array of the rectilinear grid.
106 // \param zcoords the z-coordinates array of the rectilinear grid.
107 // \param pointData pointer to the node-centered fields of the grid.
108 // \pre gridID >= 0. The code uses values of gridID < -1 as flag internally.
109 // \pre vtkStructuredExtent::Smaller(extent,wholeExtent) == true.
110 // \note A rank with no or an empty grid, should not call this method.
112 const int gridID,
113 int extent[6],
114 vtkDataArray* xcoords,
115 vtkDataArray* ycoords,
116 vtkDataArray* zcoords,
117 vtkPointData* pointData
118 );
119
127
133
142
148 void GetOutputStructuredGrid(const int gridID, vtkStructuredGrid* grid);
149
155 void GetOutputImageData(const int gridID, vtkImageData* grid);
156
162 void GetOutputRectilinearGrid(const int gridID, vtkRectilinearGrid* grid);
163
164protected:
167
169
170 vtk::detail::DomainMetaData* DomainInfo;
171 vtk::detail::StructuredGrid* InputGrid;
172 vtk::detail::StructuredGrid* OutputGrid;
173 vtk::detail::CommunicationManager* CommManager;
174
179
183 void PackData(int ext[6], vtkMultiProcessStream& bytestream);
184
188 void UnPackData(unsigned char* buffer, unsigned int size);
189
193 void AllocateBuffers(const int dim);
194
199
204
209 void GrowGrid(const int dim);
210
216 void UpdateNeighborList(const int dim);
217
222
229
230private:
232 void operator=(const vtkStructuredImplicitConnectivity&) VTK_DELETE_FUNCTION;
233};
234#endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
topologically and geometrically regular array of data
Definition: vtkImageData.h:46
a simple class to control print indentation
Definition: vtkIndent.h:40
Process communication using MPI.
stream used to pass data across processes using vtkMultiProcessController.
abstract base class for most VTK objects
Definition: vtkObject.h:60
represent and manipulate point attribute data
Definition: vtkPointData.h:38
represent and manipulate 3D points
Definition: vtkPoints.h:40
a dataset that is topologically regular with variable spacing in the three coordinate directions
topologically regular array of data
a distributed structured dataset that is implicitly connected among partitions without abutting.
void GetOutputRectilinearGrid(const int gridID, vtkRectilinearGrid *grid)
Gets the output rectilinear grid instance on this process.
void GetOutputStructuredGrid(const int gridID, vtkStructuredGrid *grid)
Gets the output structured grid instance on this process.
void ConstructOutput()
Constructs the output data-structures.
void PackData(int ext[6], vtkMultiProcessStream &bytestream)
Packs the data to send into a bytestream.
void RegisterRectilinearGrid(const int gridID, int extent[6], vtkDataArray *xcoords, vtkDataArray *ycoords, vtkDataArray *zcoords, vtkPointData *pointData)
void EstablishConnectivity()
Finds implicit connectivity for a distributed structured dataset.
static vtkStructuredImplicitConnectivity * New()
void ExchangeData()
Exchanges one layer (row or column) of data between neighboring grids to fix the implicit connectivit...
void GetOutputImageData(const int gridID, vtkImageData *grid)
Gets the output uniform grid instance on this process.
bool HasImplicitConnectivity()
Checks if there is implicit connectivity.
void SetWholeExtent(int wholeExt[6])
Sets the whole extent for the distributed structured domain.
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
void AllocateBuffers(const int dim)
Allocates send/rcv buffers needed to carry out the communication.
void ExchangeExtents()
Exchanges extents among processes.
void UnPackData(unsigned char *buffer, unsigned int size)
Unpacks the data to the output grid.
bool GlobalDataDescriptionMatch()
Checks if the data description matches globally.
void UpdateNeighborList(const int dim)
Updates the list of neighbors after growing the grid along the given dimension dim.
void GetGlobalImplicitConnectivityState()
Gets whether there is implicit connectivity across all processes.
void RegisterGrid(const int gridID, int extent[6], vtkPoints *gridPnts, vtkPointData *pointData)
void ComputeNeighbors()
Computes the neighbors with implicit connectivity.
vtk::detail::CommunicationManager * CommManager
void GrowGrid(const int dim)
Grows grid along a given dimension.
@ extent
Definition: vtkX3D.h:345
@ size
Definition: vtkX3D.h:253