VTK  9.2.6
vtkSLACReader.h
Go to the documentation of this file.
1// -*- c++ -*-
2/*=========================================================================
3
4 Program: Visualization Toolkit
5 Module: vtkSLACReader.h
6
7 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
8 All rights reserved.
9 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
10
11 This software is distributed WITHOUT ANY WARRANTY; without even
12 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
13 PURPOSE. See the above copyright notice for more information.
14
15=========================================================================*/
16
17/*-------------------------------------------------------------------------
18 Copyright 2008 Sandia Corporation.
19 Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
20 the U.S. Government retains certain rights in this software.
21-------------------------------------------------------------------------*/
22
38#ifndef vtkSLACReader_h
39#define vtkSLACReader_h
40
41#include "vtkIONetCDFModule.h" // For export macro
43
44#include "vtkSmartPointer.h" // For internal method.
45
47class vtkDoubleArray;
48class vtkIdTypeArray;
51
52class VTKIONETCDF_EXPORT vtkSLACReader : public vtkMultiBlockDataSetAlgorithm
53{
54public:
56 static vtkSLACReader* New();
57 void PrintSelf(ostream& os, vtkIndent indent) override;
58
59 vtkGetFilePathMacro(MeshFileName);
60 vtkSetFilePathMacro(MeshFileName);
61
63
68 virtual void AddModeFileName(VTK_FILEPATH const char* fname);
69 virtual void RemoveAllModeFileNames();
70 virtual unsigned int GetNumberOfModeFileNames();
71 virtual VTK_FILEPATH const char* GetModeFileName(unsigned int idx);
73
75
78 vtkGetMacro(ReadInternalVolume, vtkTypeBool);
79 vtkSetMacro(ReadInternalVolume, vtkTypeBool);
80 vtkBooleanMacro(ReadInternalVolume, vtkTypeBool);
82
84
87 vtkGetMacro(ReadExternalSurface, vtkTypeBool);
88 vtkSetMacro(ReadExternalSurface, vtkTypeBool);
89 vtkBooleanMacro(ReadExternalSurface, vtkTypeBool);
91
93
97 vtkGetMacro(ReadMidpoints, vtkTypeBool);
98 vtkSetMacro(ReadMidpoints, vtkTypeBool);
99 vtkBooleanMacro(ReadMidpoints, vtkTypeBool);
101
103
107 virtual const char* GetVariableArrayName(int index);
108 virtual int GetVariableArrayStatus(const char* name);
109 virtual void SetVariableArrayStatus(const char* name, int status);
111
113
116 virtual void ResetFrequencyScales();
117 virtual void SetFrequencyScale(int index, double scale);
119
121
124 virtual void ResetPhaseShifts();
125 virtual void SetPhaseShift(int index, double shift);
127
129
135
139 static int CanReadFile(VTK_FILEPATH const char* filename);
140
146
152
154
163
165
169 class VTKIONETCDF_EXPORT EdgeEndpoints
170 {
171 public:
173 : MinEndPoint(-1)
174 , MaxEndPoint(-1)
175 {
176 }
177 EdgeEndpoints(vtkIdType endpointA, vtkIdType endpointB)
178 {
179 if (endpointA < endpointB)
180 {
181 this->MinEndPoint = endpointA;
182 this->MaxEndPoint = endpointB;
183 }
184 else
185 {
186 this->MinEndPoint = endpointB;
187 this->MaxEndPoint = endpointA;
188 }
189 }
190 inline vtkIdType GetMinEndPoint() const { return this->MinEndPoint; }
191 inline vtkIdType GetMaxEndPoint() const { return this->MaxEndPoint; }
192 inline bool operator==(const EdgeEndpoints& other) const
193 {
194 return ((this->GetMinEndPoint() == other.GetMinEndPoint()) &&
195 (this->GetMaxEndPoint() == other.GetMaxEndPoint()));
196 }
197
198 protected:
201 };
203
205
208 class VTKIONETCDF_EXPORT MidpointCoordinates
209 {
210 public:
212 MidpointCoordinates(const double coord[3], vtkIdType id)
213 {
214 this->Coordinate[0] = coord[0];
215 this->Coordinate[1] = coord[1];
216 this->Coordinate[2] = coord[2];
217 this->ID = id;
218 }
219 double Coordinate[3];
221 };
223
224 enum
225 {
226 SURFACE_OUTPUT = 0,
227 VOLUME_OUTPUT = 1,
228 NUM_OUTPUTS = 2
229 };
230
231protected:
233 ~vtkSLACReader() override;
234
235 class vtkInternal;
236 vtkInternal* Internal;
237
238 // Friend so vtkInternal can access MidpointIdMap
239 // (so Sun CC compiler doesn't complain).
240 friend class vtkInternal;
241
243
247
252
257
262
264 vtkInformationVector* outputVector) override;
265
267 vtkInformationVector* outputVector) override;
268
273 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
274
282 virtual vtkIdType GetNumTuplesInVariable(int ncFD, int varId, int expectedNumComponents);
283
288 virtual int CheckTetrahedraWinding(int meshFD);
289
294 virtual int ReadConnectivity(
295 int meshFD, vtkMultiBlockDataSet* surfaceOutput, vtkMultiBlockDataSet* volumeOutput);
296
298
301 virtual int ReadTetrahedronInteriorArray(int meshFD, vtkIdTypeArray* connectivity);
302 virtual int ReadTetrahedronExteriorArray(int meshFD, vtkIdTypeArray* connectivity);
304
309
313 enum
314 {
315 NumPerTetInt = 5,
316 NumPerTetExt = 9
317 };
318
320
323 class VTKIONETCDF_EXPORT MidpointCoordinateMap
324 {
325 public:
329
330 void AddMidpoint(const EdgeEndpoints& edge, const MidpointCoordinates& midpoint);
331 void RemoveMidpoint(const EdgeEndpoints& edge);
334
340
341 protected:
342 class vtkInternal;
343 vtkInternal* Internal;
344
345 private:
346 // Too lazy to implement these.
348 void operator=(const MidpointCoordinateMap&) = delete;
349 };
350
352
355 class VTKIONETCDF_EXPORT MidpointIdMap
356 {
357 public:
361
362 void AddMidpoint(const EdgeEndpoints& edge, vtkIdType midpoint);
363 void RemoveMidpoint(const EdgeEndpoints& edge);
366
371
380
381 protected:
382 class vtkInternal;
383 vtkInternal* Internal;
384
385 private:
386 // Too lazy to implement these.
387 MidpointIdMap(const MidpointIdMap&) = delete;
388 void operator=(const MidpointIdMap&) = delete;
389 };
390
395 virtual int ReadCoordinates(int meshFD, vtkMultiBlockDataSet* output);
396
403 int meshFD, vtkMultiBlockDataSet* output, MidpointCoordinateMap& map);
404
410 virtual int ReadMidpointData(
411 int meshFD, vtkMultiBlockDataSet* output, MidpointIdMap& midpointIds);
412
417 virtual int RestoreMeshCache(vtkMultiBlockDataSet* surfaceOutput,
418 vtkMultiBlockDataSet* volumeOutput, vtkMultiBlockDataSet* compositeOutput);
419
424 virtual int ReadFieldData(const int* modeFDArray, int numModeFDs, vtkMultiBlockDataSet* output);
425
431
438
443 virtual int MeshUpToDate();
444
445private:
446 vtkSLACReader(const vtkSLACReader&) = delete;
447 void operator=(const vtkSLACReader&) = delete;
448};
449
450#endif // vtkSLACReader_h
Store on/off settings for data arrays, etc.
dynamic, self-adjusting array of double
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:40
Key for integer values in vtkInformation.
Key for vtkObjectBase values.
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
Composite dataset that organizes datasets into blocks.
abstract base class for most VTK objects
Definition vtkObject.h:63
Simple class used internally to define an edge based on the endpoints.
vtkIdType GetMinEndPoint() const
vtkIdType GetMaxEndPoint() const
bool operator==(const EdgeEndpoints &other) const
EdgeEndpoints(vtkIdType endpointA, vtkIdType endpointB)
Manages a map from edges to midpoint coordinates.
MidpointCoordinates * FindMidpoint(const EdgeEndpoints &edge)
Finds the coordinates for the given edge or returns nullptr if it does not exist.
void RemoveMidpoint(const EdgeEndpoints &edge)
void AddMidpoint(const EdgeEndpoints &edge, const MidpointCoordinates &midpoint)
Simple class used internally for holding midpoint information.
MidpointCoordinates(const double coord[3], vtkIdType id)
Manages a map from edges to the point id of the midpoint.
vtkIdType GetNumberOfMidpoints() const
void AddMidpoint(const EdgeEndpoints &edge, vtkIdType midpoint)
void InitTraversal()
Initialize iteration.
bool GetNextMidpoint(EdgeEndpoints &edge, vtkIdType &midpoint)
Get the next midpoint in the iteration.
vtkIdType * FindMidpoint(const EdgeEndpoints &edge)
Finds the id for the given edge or returns nullptr if it does not exist.
void RemoveMidpoint(const EdgeEndpoints &edge)
A reader for a data format used by Omega3p, Tau3p, and several other tools used at the Standford Line...
virtual void ResetFrequencyScales()
Sets the scale factor for each mode.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
bool TimeStepModes
True if "mode" files are a sequence of time steps.
vtkInternal * Internal
virtual void SetPhaseShift(int index, double shift)
Sets the phase offset for each mode.
virtual void SetFrequencyScale(int index, double scale)
Sets the scale factor for each mode.
virtual int ReadMidpointCoordinates(int meshFD, vtkMultiBlockDataSet *output, MidpointCoordinateMap &map)
Reads in the midpoint coordinate data from the mesh file and returns a map from edges to midpoints.
virtual const char * GetVariableArrayName(int index)
Variable array selection.
static vtkInformationIntegerKey * IS_EXTERNAL_SURFACE()
This key is attached to the metadata information of all data sets in the output that are part of the ...
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Callback registered with the VariableArraySelection.
static int CanReadFile(VTK_FILEPATH const char *filename)
Returns true if the given file can be read by this reader.
static vtkSLACReader * New()
vtkTypeBool ReadMidpoints
static vtkInformationObjectBaseKey * POINT_DATA()
All the data sets stored in the multiblock output share the same point data.
virtual VTK_FILEPATH const char * GetModeFileName(unsigned int idx)
There may be one mode file (usually for actual modes) or multiple mode files (which usually actually ...
~vtkSLACReader() override
vtkGetFilePathMacro(MeshFileName)
virtual int RestoreMeshCache(vtkMultiBlockDataSet *surfaceOutput, vtkMultiBlockDataSet *volumeOutput, vtkMultiBlockDataSet *compositeOutput)
Instead of reading data from the mesh file, restore the data from the previous mesh file read.
virtual int ReadTetrahedronExteriorArray(int meshFD, vtkIdTypeArray *connectivity)
Reads tetrahedron connectivity arrays.
virtual int ReadFieldData(const int *modeFDArray, int numModeFDs, vtkMultiBlockDataSet *output)
Read in the field data from the mode file.
virtual void ResetPhaseShifts()
Sets the phase offset for each mode.
virtual int MeshUpToDate()
Returns 1 if the mesh is up to date, 0 if the mesh needs to be read from disk.
virtual vtkDoubleArray * GetPhaseShifts()
NOTE: This is not thread-safe.
virtual int ReadMidpointData(int meshFD, vtkMultiBlockDataSet *output, MidpointIdMap &midpointIds)
Read in the midpoint data from the mesh file.
virtual unsigned int GetNumberOfModeFileNames()
There may be one mode file (usually for actual modes) or multiple mode files (which usually actually ...
bool FrequencyModes
True if mode files describe vibrating fields.
virtual int ReadCoordinates(int meshFD, vtkMultiBlockDataSet *output)
Read in the point coordinate data from the mesh file.
virtual vtkDoubleArray * GetFrequencyScales()
NOTE: This is not thread-safe.
virtual int ReadTetrahedronInteriorArray(int meshFD, vtkIdTypeArray *connectivity)
Reads tetrahedron connectivity arrays.
vtkTimeStamp MeshReadTime
A time stamp for the last time the mesh file was read.
virtual int ReadConnectivity(int meshFD, vtkMultiBlockDataSet *surfaceOutput, vtkMultiBlockDataSet *volumeOutput)
Read the connectivity information from the mesh file.
bool ReadModeData
True if reading from a proper mode file.
static vtkInformationObjectBaseKey * POINTS()
All the data sets stored in the multiblock output share the same point data.
virtual int GetVariableArrayStatus(const char *name)
Variable array selection.
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
virtual vtkIdType GetNumTuplesInVariable(int ncFD, int varId, int expectedNumComponents)
Convenience function that checks the dimensions of a 2D netCDF array that is supposed to be a set of ...
vtkTypeBool ReadInternalVolume
vtkTypeBool ReadExternalSurface
static vtkInformationIntegerKey * IS_INTERNAL_VOLUME()
This key is attached to the metadata information of all data sets in the output that are part of the ...
virtual int InterpolateMidpointData(vtkMultiBlockDataSet *output, MidpointIdMap &map)
Takes the data read on the fields and interpolates data for the midpoints.
virtual void AddModeFileName(VTK_FILEPATH const char *fname)
There may be one mode file (usually for actual modes) or multiple mode files (which usually actually ...
virtual int GetNumberOfVariableArrays()
Variable array selection.
virtual void RemoveAllModeFileNames()
There may be one mode file (usually for actual modes) or multiple mode files (which usually actually ...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSetFilePathMacro(MeshFileName)
virtual int CheckTetrahedraWinding(int meshFD)
Checks the winding of the tetrahedra in the mesh file.
virtual void SetVariableArrayStatus(const char *name, int status)
Variable array selection.
virtual vtkSmartPointer< vtkDataArray > ReadPointDataArray(int ncFD, int varId)
Reads point data arrays.
Hold a reference to a vtkObjectBase instance.
record modification and/or execution time
int vtkTypeBool
Definition vtkABI.h:69
int vtkIdType
Definition vtkType.h:332
#define VTK_FILEPATH