1 /*=========================================================================
2 
3   Program:   ParaView
4   Module:    vtkModelMetadata.h
5 
6   Copyright (c) Kitware, Inc.
7   All rights reserved.
8   See Copyright.txt or http://www.paraview.org/HTML/Copyright.html 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 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
20 /**
21  * @class   vtkModelMetadata
22  * @brief   This class encapsulates the metadata
23  *   that appear in mesh-based file formats but do not appear in
24  *   vtkUnstructuredGrid.
25  *
26  *
27  *   This class is inspired by the Exodus II file format, but
28  *   because this class does not depend on the Exodus library, it
29  *   should be possible to use it to represent metadata for other
30  *   dataset file formats.  Sandia Labs uses it in their Exodus II
31  *   reader, their Exodus II writer and their EnSight writer.
32  *   vtkDistributedDataFilter looks for metadata attached to
33  *   it's input and redistributes the metadata with the grid.
34  *
35  *   The fields in this class are those described in the document
36  *   "EXODUS II: A Finite Element Data Model", SAND92-2137, November 1995.
37  *
38  *   Element and node IDs stored in this object must be global IDs,
39  *   in the event that the original dataset was partitioned across
40  *   many files.
41  *
42  *   One way to initialize this object is by using vtkExodusModel
43  *   (a Sandia class used by the Sandia Exodus reader).
44  *   That class will take an open Exodus II file and a
45  *   vtkUnstructuredGrid drawn from it and will set the required fields.
46  *
47  *   Alternatively, you can use all the Set*
48  *   methods to set the individual fields. This class does not
49  *   copy the data, it simply uses your pointer. This
50  *   class will free the storage associated with your pointer
51  *   when the class is deleted.  Most fields have sensible defaults.
52  *   The only requirement is that if you are using this ModelMetadata
53  *   to write out an Exodus or EnSight file in parallel, you must
54  *   SetBlockIds and SetBlockIdArrayName.  Your vtkUnstructuredGrid must
55  *   have a cell array giving the block ID for each cell.
56  *
57  * @warning
58  *   The Exodus II library supports an optimized element order map
59  *   (section 3.7 in the SAND document).  It contains all the element
60  *   IDs, listed in the order in which a solver should process them.
61  *   We don't include this, and won't unless there is a request.
62  *
63  * @warning
64  *   There is an assumption in some classes that the name of the cell
65  *   array containing global element ids is "GlobalElementId" and the
66  *   name of the point array containing global node ids is "GlobalNodeId".
67  *   (element == cell) and (node == point).
68  *
69  * @sa
70  *   vtkDistributedDataFilter vtkExtractCells
71  */
72 
73 #ifndef vtkModelMetadata_h
74 #define vtkModelMetadata_h
75 
76 #include "vtkIOExodusModule.h" // For export macro
77 #include "vtkObject.h"
78 #include "vtkSmartPointer.h" // for vtkSmartPointer
79 #include "vtkStringArray.h"  // for vtkStringArray
80 class vtkDataSet;
81 class vtkCharArray;
82 class vtkIdTypeArray;
83 class vtkIntArray;
84 class vtkFloatArray;
85 class vtkIntArray;
86 class vtkStringArray;
87 class vtkModelMetadataSTLCloak;
88 
89 class VTKIOEXODUS_EXPORT vtkModelMetadata : public vtkObject
90 {
91 public:
92   vtkTypeMacro(vtkModelMetadata, vtkObject);
93   void PrintSelf(ostream& os, vtkIndent indent) override;
94   static vtkModelMetadata* New();
95 
96   /**
97    * The global fields are those which pertain to the whole
98    * file.  Examples are the title, information lines,
99    * and list of block IDs.  This method prints out all the
100    * global information.
101    */
102 
103   virtual void PrintGlobalInformation();
104 
105   /**
106    * The local fields are those which depend on exactly which
107    * blocks, which time step, and which variables you read in
108    * from the file.  Examples are the number of cells in
109    * each block, and the list of nodes in a node set, or the
110    * value of the global variables at a time step.  If
111    * VERBOSE_TESTING is defined in your execution environment,
112    * this method will print more than mere counts, and actually
113    * print a few of the IDs, distribution factors and so on.  If
114    * VERY_VERBOSE_TESTING is defined, it will print out
115    * all ID lists, distribution factor lists, and so on.
116    */
117 
118   virtual void PrintLocalInformation();
119 
120   ///@{
121   /**
122    * The title of the dataset.
123    */
124   vtkSetStringMacro(Title);
GetTitle()125   const char* GetTitle() const { return this->Title; }
126   ///@}
127 
128   /**
129    * Set the information lines.
130    */
131   void SetInformationLines(int numLines, char** lines);
132 
133   /**
134    * Get a pointer to all the information lines.  The number
135    * of lines is returned;
136    */
137   int GetInformationLines(char*** lines) const;
138 
139   /**
140    * Get the number of information lines.
141    */
GetNumberOfInformationLines()142   int GetNumberOfInformationLines() const { return this->NumberOfInformationLines; }
143 
144   ///@{
145   /**
146    * Set the index of the time step represented by the results
147    * data in the file attached to this ModelMetadata object.  Time
148    * step indices start at 0 in this file, they start at 1 in
149    * an Exodus file.
150    */
151   vtkSetMacro(TimeStepIndex, int);
GetTimeStepIndex()152   int GetTimeStepIndex() const { return this->TimeStepIndex; }
153   ///@}
154 
155   /**
156    * Set the total number of time steps in the file,
157    * and the value at each time step.  We use your time
158    * step value array and delete it when we're done.
159    */
160   void SetTimeSteps(int numberOfTimeSteps, float* timeStepValues);
GetNumberOfTimeSteps()161   int GetNumberOfTimeSteps() const { return this->NumberOfTimeSteps; }
162 
163   /**
164    * Get the time step values
165    */
GetTimeStepValues()166   float* GetTimeStepValues() const { return this->TimeStepValues; }
167 
168   /**
169    * The name of the one, two or three coordinate dimensions.
170    */
171   void SetCoordinateNames(int dimension, char**);
GetCoordinateNames()172   char** GetCoordinateNames() const { return this->CoordinateNames; }
173 
174   /**
175    * Get the dimension of the model.  This is also the number
176    * of coordinate names.
177    */
GetDimension()178   int GetDimension() const { return this->Dimension; }
179 
180   ///@{
181   /**
182    * The number of blocks in the file.  Set this before setting
183    * any of the block arrays.
184    */
185   vtkSetMacro(NumberOfBlocks, int);
GetNumberOfBlocks()186   int GetNumberOfBlocks() const { return this->NumberOfBlocks; }
187   ///@}
188 
189   /**
190    * An arbitrary integer ID for each block.
191    * We use your pointer, and free the memory when the object is freed.
192    */
193   void SetBlockIds(int*);
GetBlockIds()194   int* GetBlockIds() const { return this->BlockIds; }
195 
196   /**
197    * Element type for each block - a name that means
198    * something to person who created the file.
199    * We use your pointers, and free the memory when the object is freed.
200    */
201   void SetBlockElementType(char**);
GetBlockElementType()202   char** GetBlockElementType() const { return this->BlockElementType; }
203 
204   /**
205    * Set or get a pointer to a list of the number of elements in
206    * each block.
207    * We use your pointers, and free the memory when the object is freed.
208    */
209   int SetBlockNumberOfElements(int* nelts);
GetBlockNumberOfElements()210   int* GetBlockNumberOfElements() const { return this->BlockNumberOfElements; }
211 
212   /**
213    * Set or get a pointer to a list of the number of nodes in the
214    * elements of each block.
215    * We use your pointers, and free the memory when the object is freed.
216    */
217   void SetBlockNodesPerElement(int*);
GetBlockNodesPerElement()218   int* GetBlockNodesPerElement() const { return this->BlockNodesPerElement; }
219 
220   /**
221    * Set or get a pointer to a list global element IDs for the
222    * elements in each block.
223    * We use your pointers, and free the memory when the object is freed.
224    */
225   void SetBlockElementIdList(int*);
GetBlockElementIdList()226   int* GetBlockElementIdList() const { return this->BlockElementIdList; }
227 
228   /**
229    * Get the length of the list of elements in every block.
230    */
GetSumElementsPerBlock()231   int GetSumElementsPerBlock() const { return this->SumElementsPerBlock; }
232 
233   /**
234    * Get a list of the index into the BlockElementIdList of the
235    * start of each block's elements.
236    */
GetBlockElementIdListIndex()237   int* GetBlockElementIdListIndex() const { return this->BlockElementIdListIndex; }
238 
239   /**
240    * Set or get a pointer to a list of the number of attributes
241    * stored for the elements in each block.
242    * We use your pointers, and free the memory when the object is freed.
243    */
244   int SetBlockNumberOfAttributesPerElement(int* natts);
GetBlockNumberOfAttributesPerElement()245   int* GetBlockNumberOfAttributesPerElement() const
246   {
247     return this->BlockNumberOfAttributesPerElement;
248   }
249 
250   /**
251    * Set or get a pointer to a list of the attributes for all
252    * blocks.  The order of the list should be by block, by element
253    * within the block, by attribute.  Omit blocks that don't
254    * have element attributes.
255    */
256   void SetBlockAttributes(float*);
GetBlockAttributes()257   float* GetBlockAttributes() const { return this->BlockAttributes; }
258 
259   /**
260    * Get the length of the list of floating point block attributes.
261    */
GetSizeBlockAttributeArray()262   int GetSizeBlockAttributeArray() const { return this->SizeBlockAttributeArray; }
263 
264   /**
265    * Get a list of the index into the BlockAttributes of the
266    * start of each block's element attribute list.
267    */
GetBlockAttributesIndex()268   int* GetBlockAttributesIndex() const { return this->BlockAttributesIndex; }
269 
270   ///@{
271   /**
272    * The number of node sets in the file.  Set this value before
273    * setting the various node set arrays.
274    */
275   vtkSetMacro(NumberOfNodeSets, int);
GetNumberOfNodeSets()276   int GetNumberOfNodeSets() const { return this->NumberOfNodeSets; }
277   ///@}
278 
SetNodeSetNames(vtkStringArray * names)279   void SetNodeSetNames(vtkStringArray* names) { this->NodeSetNames = names; }
GetNodeSetNames()280   vtkStringArray* GetNodeSetNames() const { return this->NodeSetNames; }
281 
282   /**
283    * Set or get the list the IDs for each node set.
284    * Length of list is the number of node sets.
285    * We use your pointer, and free the memory when the object is freed.
286    */
287   void SetNodeSetIds(int*);
GetNodeSetIds()288   int* GetNodeSetIds() const { return this->NodeSetIds; }
289 
290   /**
291    * Set or get a pointer to a list of the number of nodes in each node set.
292    * We use your pointer, and free the memory when the object is freed.
293    */
294   void SetNodeSetSize(int*);
GetNodeSetSize()295   int* GetNodeSetSize() const { return this->NodeSetSize; }
296 
297   /**
298    * Set or get a pointer to a concatenated list of the
299    * IDs of all nodes in each node set.  First list all IDs in
300    * node set 0, then all IDs in node set 1, and so on.
301    * We use your pointer, and free the memory when the object is freed.
302    */
303   void SetNodeSetNodeIdList(int*);
GetNodeSetNodeIdList()304   int* GetNodeSetNodeIdList() const { return this->NodeSetNodeIdList; }
305 
306   /**
307    * Set or get a list of the number of distribution factors stored
308    * by each node set.  This is either 0 or equal to the number of
309    * nodes in the node set.
310    * Length of list is number of node sets.
311    * We use your pointer, and free the memory when the object is freed.
312    */
313   void SetNodeSetNumberOfDistributionFactors(int*);
GetNodeSetNumberOfDistributionFactors()314   int* GetNodeSetNumberOfDistributionFactors() const
315   {
316     return this->NodeSetNumberOfDistributionFactors;
317   }
318 
319   /**
320    * Set or get a list of the distribution factors for the node sets.
321    * The list is organized by node set, and within node set by node.
322    * We use your pointer, and free the memory when the object is freed.
323    */
324   void SetNodeSetDistributionFactors(float*);
GetNodeSetDistributionFactors()325   float* GetNodeSetDistributionFactors() const { return this->NodeSetDistributionFactors; }
326 
327   ///@{
328   /**
329    * Get the total number of nodes in all node sets
330    */
331   vtkSetMacro(SumNodesPerNodeSet, int);
GetSumNodesPerNodeSet()332   int GetSumNodesPerNodeSet() const { return this->SumNodesPerNodeSet; }
333   ///@}
334 
335   /**
336    * Get the total number of distribution factors stored for all node sets
337    */
GetSumDistFactPerNodeSet()338   int GetSumDistFactPerNodeSet() const { return this->SumDistFactPerNodeSet; }
339 
340   /**
341    * Get a list of the index of the starting entry for each node set
342    * in the list of node set node IDs.
343    */
GetNodeSetNodeIdListIndex()344   int* GetNodeSetNodeIdListIndex() const { return this->NodeSetNodeIdListIndex; }
345 
346   /**
347    * Get a list of the index of the starting entry for each node set
348    * in the list of node set distribution factors.
349    */
GetNodeSetDistributionFactorIndex()350   int* GetNodeSetDistributionFactorIndex() const { return this->NodeSetDistributionFactorIndex; }
351 
352   ///@{
353   /**
354    * Set or get the number of side sets.  Set this value before
355    * setting any of the other side set arrays.
356    */
357   vtkSetMacro(NumberOfSideSets, int);
GetNumberOfSideSets()358   int GetNumberOfSideSets() const { return this->NumberOfSideSets; }
359   ///@}
360 
SetSideSetNames(vtkStringArray * names)361   void SetSideSetNames(vtkStringArray* names) { this->SideSetNames = names; }
GetSideSetNames()362   vtkStringArray* GetSideSetNames() const { return this->SideSetNames; }
363 
364   /**
365    * Set or get a pointer to a list giving the ID of each side set.
366    * We use your pointer, and free the memory when the object is freed.
367    */
368   void SetSideSetIds(int*);
GetSideSetIds()369   int* GetSideSetIds() const { return this->SideSetIds; }
370 
371   /**
372    * Set or get a pointer to a list of the number of sides in each side set.
373    * We use your pointer, and free the memory when the object is freed.
374    */
375   int SetSideSetSize(int* sizes);
GetSideSetSize()376   int* GetSideSetSize() const { return this->SideSetSize; }
377 
378   /**
379    * Set or get a pointer to a list of the number of distribution
380    * factors stored by each side set.   Each side set has either
381    * no distribution factors, or 1 per node in the side set.
382    * We use your pointer, and free the memory when the object is freed.
383    */
384   int SetSideSetNumberOfDistributionFactors(int* df);
GetSideSetNumberOfDistributionFactors()385   int* GetSideSetNumberOfDistributionFactors() const
386   {
387     return this->SideSetNumberOfDistributionFactors;
388   }
389 
390   /**
391    * Set or get a pointer to a list of the elements containing each
392    * side in each side set.  The list is organized by side set, and
393    * within side set by element.
394    * We use your pointer, and free the memory when the object is freed.
395    */
396   void SetSideSetElementList(int*);
GetSideSetElementList()397   int* GetSideSetElementList() const { return this->SideSetElementList; }
398 
399   /**
400    * Set or get a pointer to the element side for each side in the side set.
401    * (See the manual for the convention for numbering sides in different
402    * types of cells.)  Side Ids are arranged by side set and within
403    * side set by side, and correspond to the SideSetElementList.
404    * We use your pointer, and free the memory when the object is freed.
405    */
406   void SetSideSetSideList(int*);
GetSideSetSideList()407   int* GetSideSetSideList() const { return this->SideSetSideList; }
408 
409   /**
410    * Set or get a pointer to a list of the number of nodes in each
411    * side of each side set.  This list is organized by side set, and
412    * within side set by side.
413    * We use your pointer, and free the memory when the object is freed.
414    */
415   void SetSideSetNumDFPerSide(int* numNodes);
GetSideSetNumDFPerSide()416   int* GetSideSetNumDFPerSide() const { return this->SideSetNumDFPerSide; }
417 
418   /**
419    * Set or get a pointer to a list of all the distribution factors.
420    * For every side set that has distribution factors, the number of
421    * factors per node was given in the SideSetNumberOfDistributionFactors
422    * array.  If this number for a given side set is N, then for that
423    * side set we have N floating point values for each node for each
424    * side in the side set.  If nodes are repeated in more than one
425    * side, we repeat the distribution factors.  So this list is in order
426    * by side set, by node.
427    * We use your pointer, and free the memory when the object is freed.
428    */
429   void SetSideSetDistributionFactors(float*);
GetSideSetDistributionFactors()430   float* GetSideSetDistributionFactors() const { return this->SideSetDistributionFactors; }
431 
432   ///@{
433   /**
434    * Get the total number of sides in all side sets
435    */
436   vtkSetMacro(SumSidesPerSideSet, int);
GetSumSidesPerSideSet()437   int GetSumSidesPerSideSet() const { return this->SumSidesPerSideSet; }
438   ///@}
439 
440   /**
441    * Get the total number of distribution factors stored for all side sets
442    */
GetSumDistFactPerSideSet()443   int GetSumDistFactPerSideSet() const { return this->SumDistFactPerSideSet; }
444 
445   /**
446    * Get a list of the index of the starting entry for each side set
447    * in the list of side set side IDs.
448    */
GetSideSetListIndex()449   int* GetSideSetListIndex() const { return this->SideSetListIndex; }
450 
451   /**
452    * Get a list of the index of the starting entry for each side set
453    * in the list of side set distribution factors.
454    */
GetSideSetDistributionFactorIndex()455   int* GetSideSetDistributionFactorIndex() const { return this->SideSetDistributionFactorIndex; }
456 
457   /**
458    * The number of block properties (global variables)
459    */
GetNumberOfBlockProperties()460   int GetNumberOfBlockProperties() const { return this->NumberOfBlockProperties; }
461 
462   /**
463    * Set or get the names of the block properties.
464    */
465   void SetBlockPropertyNames(int numProp, char** names);
GetBlockPropertyNames()466   char** GetBlockPropertyNames() const { return this->BlockPropertyNames; }
467 
468   /**
469    * Set or get value for each variable for each block.  List
470    * the integer values in order by variable and within variable
471    * by block.
472    */
473   void SetBlockPropertyValue(int*);
GetBlockPropertyValue()474   int* GetBlockPropertyValue() const { return this->BlockPropertyValue; }
475 
476   /**
477    * The number of node set properties (global variables)
478    */
GetNumberOfNodeSetProperties()479   int GetNumberOfNodeSetProperties() const { return this->NumberOfNodeSetProperties; }
480 
481   /**
482    * Set or get the names of the node setproperties.
483    */
484   void SetNodeSetPropertyNames(int numProp, char** names);
GetNodeSetPropertyNames()485   char** GetNodeSetPropertyNames() const { return this->NodeSetPropertyNames; }
486 
487   /**
488    * Set or get value for each variable for each node set.  List
489    * the integer values in order by variable and within variable
490    * by node set.
491    */
492   void SetNodeSetPropertyValue(int*);
GetNodeSetPropertyValue()493   int* GetNodeSetPropertyValue() const { return this->NodeSetPropertyValue; }
494 
495   /**
496    * The number of side set properties (global variables)
497    */
GetNumberOfSideSetProperties()498   int GetNumberOfSideSetProperties() const { return this->NumberOfSideSetProperties; }
499 
500   /**
501    * Set or get the names of the side set properties.
502    */
503   void SetSideSetPropertyNames(int numProp, char** names);
GetSideSetPropertyNames()504   char** GetSideSetPropertyNames() const { return this->SideSetPropertyNames; }
505 
506   /**
507    * Set or get value for each variable for each side set.  List
508    * the integer values in order by variable and within variable
509    * by side set.
510    */
511   void SetSideSetPropertyValue(int*);
GetSideSetPropertyValue()512   int* GetSideSetPropertyValue() const { return this->SideSetPropertyValue; }
513 
514   /**
515    * Get the number of global variables per time step
516    */
GetNumberOfGlobalVariables()517   int GetNumberOfGlobalVariables() const { return this->NumberOfGlobalVariables; }
518 
519   /**
520    * Set or get the names of the global variables
521    */
522   void SetGlobalVariableNames(int numVarNames, char** n);
GetGlobalVariableNames()523   char** GetGlobalVariableNames() const { return this->GlobalVariableNames; }
524 
525   /**
526    * Set or get the values of the global variables at the current
527    * time step.
528    */
529   void SetGlobalVariableValue(float* f);
GetGlobalVariableValue()530   float* GetGlobalVariableValue() const { return this->GlobalVariableValue; }
531 
532   /**
533    * The ModelMetadata maintains a list of the element variables that
534    * were in the original file, and a list of the cell variables
535    * in the UGrid derived from that file.  Some of the scalar variables
536    * in the original file were combined into vectors in the UGrid.
537    * In this method, provide the number of original element variables,
538    * the names of the original element variables, the number of
539    * element variables in the UGrid, the number of components for each
540    * of those variables, and a map from each UGrid variable to the
541    * the variable in the list of original names that represents it's
542    * first component.
543    */
544   void SetElementVariableInfo(
545     int numOrigNames, char** origNames, int numNames, char** names, int* numComp, int* map);
546 
547   /**
548    * The ModelMetadata maintains a list of the node variables that
549    * were in the original file, and a list of the node variables
550    * in the UGrid derived from that file.  Some of the scalar variables
551    * in the original file were combined into vectors in the UGrid.
552    * In this method, provide the number of original node variables,
553    * the names of the original node variables, the number of
554    * node variables in the UGrid, the number of components for each
555    * of those variables, and a map from each UGrid variable to the
556    * the variable in the list of original names that represents it's
557    * first component.
558    */
559   void SetNodeVariableInfo(
560     int numOrigNames, char** origNames, int numNames, char** names, int* numComp, int* map);
561 
562   /**
563    * A truth table indicating which element variables are
564    * defined for which blocks. The variables are all the original
565    * element variables that were in the file.
566    * The table is by block ID and within block ID by variable.
567    */
568   void SetElementVariableTruthTable(int*);
GetElementVariableTruthTable()569   int* GetElementVariableTruthTable() const { return this->ElementVariableTruthTable; }
570 
571   ///@{
572   /**
573    * Instead of a truth table of all "1"s, you can set this
574    * instance variable to indicate that all variables are
575    * defined in all blocks.
576    */
577   vtkSetMacro(AllVariablesDefinedInAllBlocks, vtkTypeBool);
578   vtkBooleanMacro(AllVariablesDefinedInAllBlocks, vtkTypeBool);
GetAllVariablesDefinedInAllBlocks()579   vtkTypeBool GetAllVariablesDefinedInAllBlocks() const
580   {
581     return this->AllVariablesDefinedInAllBlocks;
582   }
583   ///@}
584 
585   /**
586    * The ModelMetadata object may contain these lists:
587    * o  the variables in the original data file
588    * o  the variables created in the u grid from those original variables
589    * o  a mapping from the grid variable names to the original names
590    * o  a list of the number of components each grid variable has
591 
592    * (Example: Variables in Exodus II files are all scalars.  Some are
593    * combined by the ExodusReader into vector variables in the grid.)
594 
595    * These methods return names of the original variables, the names
596    * of the grid variables, a list of the number of components in
597    * each grid variable, and a list of the index into the list of
598    * original variable names where the original name of the first
599    * component of a grid variable may be found.  The names of subsequent
600    * components would immediately follow the name of the first
601    * component.
602    */
GetOriginalNumberOfElementVariables()603   int GetOriginalNumberOfElementVariables() const { return this->OriginalNumberOfElementVariables; }
GetOriginalElementVariableNames()604   char** GetOriginalElementVariableNames() const { return this->OriginalElementVariableNames; }
GetNumberOfElementVariables()605   int GetNumberOfElementVariables() const { return this->NumberOfElementVariables; }
GetElementVariableNames()606   char** GetElementVariableNames() const { return this->ElementVariableNames; }
GetElementVariableNumberOfComponents()607   int* GetElementVariableNumberOfComponents() const
608   {
609     return this->ElementVariableNumberOfComponents;
610   }
GetMapToOriginalElementVariableNames()611   int* GetMapToOriginalElementVariableNames() const
612   {
613     return this->MapToOriginalElementVariableNames;
614   }
615 
GetOriginalNumberOfNodeVariables()616   int GetOriginalNumberOfNodeVariables() const { return this->OriginalNumberOfNodeVariables; }
GetOriginalNodeVariableNames()617   char** GetOriginalNodeVariableNames() const { return this->OriginalNodeVariableNames; }
GetNumberOfNodeVariables()618   int GetNumberOfNodeVariables() const { return this->NumberOfNodeVariables; }
GetNodeVariableNames()619   char** GetNodeVariableNames() const { return this->NodeVariableNames; }
GetNodeVariableNumberOfComponents()620   int* GetNodeVariableNumberOfComponents() const { return this->NodeVariableNumberOfComponents; }
GetMapToOriginalNodeVariableNames()621   int* GetMapToOriginalNodeVariableNames() const { return this->MapToOriginalNodeVariableNames; }
622 
623   ///@{
624   /**
625    * Free selected portions of the metadata when updating values
626    * in the vtkModelMetadata object.  Resetting a particular field,
627    * (i.e. SetNodeSetIds) frees the previous setting, but if you
628    * are not setting every field, you may want to do a wholesale
629    * "Free" first.
630 
631    * FreeAllGlobalData frees all the fields which don't depend on
632    * which time step, which blocks, or which variables are in the input.
633    * FreeAllLocalData frees all the fields which do depend on which
634    * time step, blocks or variables are in the input.
635    * FreeBlockDependentData frees all metadata fields which depend on
636    * which blocks were read in.
637    */
638   void FreeAllGlobalData();
639   void FreeAllLocalData();
640   void FreeBlockDependentData();
641   void FreeOriginalElementVariableNames();
642   void FreeOriginalNodeVariableNames();
643   void FreeUsedElementVariableNames();
644   void FreeUsedNodeVariableNames();
645   void FreeUsedElementVariables();
646   void FreeUsedNodeVariables();
647   ///@}
648 
649   /**
650    * Set the object back to it's initial state
651    */
652   void Reset();
653 
654 protected:
655   vtkModelMetadata();
656   ~vtkModelMetadata() override;
657 
658 private:
659   void InitializeAllMetadata();
660   void InitializeAllIvars();
661 
662   void FreeAllMetadata();
663   void FreeAllIvars();
664 
665   int BuildBlockElementIdListIndex();
666   int BuildBlockAttributesIndex();
667   int BuildSideSetDistributionFactorIndex();
668 
669   static char* StrDupWithNew(const char* s);
670 
671   static int FindNameOnList(char* name, char** list, int listLen);
672 
673   void ShowFloats(const char* what, int num, float* f);
674   void ShowLines(const char* what, int num, char** l);
675   void ShowIntArray(const char* what, int numx, int numy, int* id);
676   void ShowInts(const char* what, int num, int* id);
677   void ShowListsOfInts(const char* what, int* list, int nlists, int* idx, int len, int verbose);
678   void ShowListsOfFloats(const char* what, float* list, int nlists, int* idx, int len, int verbose);
679 
680   void SetOriginalElementVariableNames(int nvars, char** names);
681   void SetElementVariableNames(int nvars, char** names);
682   void SetElementVariableNumberOfComponents(int* comp);
683   void SetMapToOriginalElementVariableNames(int* map);
684 
685   void SetOriginalNodeVariableNames(int nvars, char** names);
686   void SetNodeVariableNames(int nvars, char** names);
687   void SetNodeVariableNumberOfComponents(int* comp);
688   void SetMapToOriginalNodeVariableNames(int* map);
689 
690   int CalculateMaximumLengths(int& maxString, int& maxLine);
691 
692   // Fields in Exodus II file and their size (defined in exodusII.h)
693   //   (G - global fields, relevant to entire file or file set)
694   //   (L - local fields, they differ depending on which cells and nodes are
695   //        in a file of a partitioned set, or are read in from file)
696 
697   char* Title; // (G)
698 
699   int NumberOfInformationLines; // (G)
700   char** InformationLine;       // (G)
701 
702   int Dimension;          // (G)
703   char** CoordinateNames; // (at most 3 of these) (G)
704 
705   // Time steps
706 
707   int TimeStepIndex;     // starting at 0 (Exodus file starts at 1)
708   int NumberOfTimeSteps; // (G)
709   float* TimeStepValues; // (G)
710 
711   // Block information - arrays that are input with Set*
712 
713   int NumberOfBlocks; // (G)
714 
715   int* BlockIds;                          // NumberOfBlocks (G) (start at 1)
716   char** BlockElementType;                // NumberOfBlocks (G)
717   int* BlockNumberOfElements;             // NumberOfBlocks (L)
718   int* BlockNodesPerElement;              // NumberOfBlocks (G)
719   int* BlockNumberOfAttributesPerElement; // NumberOfBlocks (G)
720   int* BlockElementIdList;                // SumElementsPerBlock     (L)
721   float* BlockAttributes;                 // SizeBlockAttributeArray (L)
722 
723   // Block information - values that we calculate
724 
725   int SumElementsPerBlock;
726   int SizeBlockAttributeArray;
727 
728   int* BlockElementIdListIndex; // NumberOfBlocks
729   int* BlockAttributesIndex;    // NumberOfBlocks
730 
731   vtkModelMetadataSTLCloak* BlockIdIndex; // computed map
732 
733   // Node Sets - arrays that are input to the class with Set*
734 
735   int NumberOfNodeSets; // (G)
736 
737   vtkSmartPointer<vtkStringArray> NodeSetNames;
738 
739   int* NodeSetIds;                         // NumberOfNodeSets (G)
740   int* NodeSetSize;                        // NumberOfNodeSets (L)
741   int* NodeSetNumberOfDistributionFactors; // NNS (L) (NSNDF[i] is 0 or NSS[i])
742   int* NodeSetNodeIdList;                  // SumNodesPerNodeSet (L)
743   float* NodeSetDistributionFactors;       // SumDistFactPerNodeSet (L)
744 
745   // Node Sets - values or arrays that the class computes
746 
747   int SumNodesPerNodeSet;
748   int SumDistFactPerNodeSet;
749 
750   int* NodeSetNodeIdListIndex;         // NumberOfNodeSets
751   int* NodeSetDistributionFactorIndex; // NumberOfNodeSets
752 
753   // Side Sets - input to class with Set*
754 
755   int NumberOfSideSets; // (G)
756 
757   vtkSmartPointer<vtkStringArray> SideSetNames;
758 
759   int* SideSetIds;                         // NumberOfSideSets (G)
760   int* SideSetSize;                        // NumberOfSideSets (L)
761   int* SideSetNumberOfDistributionFactors; // NSS (L) (SSNDF[i] = 0 or NumNodesInSide)
762   int* SideSetElementList;                 // SumSidesPerSideSet (L)
763   int* SideSetSideList;                    // SumSidesPerSideSet (L)
764   int* SideSetNumDFPerSide;                // SumSidesPerSideSet (L)
765   float* SideSetDistributionFactors;       // SumDistFactPerSideSet (L)
766 
767   // Side Sets - calculated by class
768 
769   int SumSidesPerSideSet;
770   int SumDistFactPerSideSet;
771 
772   int* SideSetListIndex;               // NumberOfSideSets
773   int* SideSetDistributionFactorIndex; // NumberOfSideSets
774 
775   // Other properties, provided as input with Set*
776 
777   int NumberOfBlockProperties; // (G)
778   char** BlockPropertyNames;   // one per property (G)
779   int* BlockPropertyValue;     // NumBlocks * NumBlockProperties (G)
780 
781   int NumberOfNodeSetProperties; // (G)
782   char** NodeSetPropertyNames;   // one per property (G)
783   int* NodeSetPropertyValue;     // NumNodeSets * NumNodeSetProperties (G)
784 
785   int NumberOfSideSetProperties; // (G)
786   char** SideSetPropertyNames;   // one per property (G)
787   int* SideSetPropertyValue;     // NumSideSets * NumSideSetProperties (G)
788 
789   // Global variables, 1 value per time step per variable.  We store
790   // these as floats, even if they are doubles in the file.  The values
791   // are global in the sense that they apply to the whole data set, but
792   // the are local in the sense that they can change with each time step.
793   // For the purpose of this object, which represents a particular
794   // time step, they are therefore considered "local".  (Since they need
795   // to be updated every time another read is done from the file.)
796 
797   int NumberOfGlobalVariables; // (G)
798   char** GlobalVariableNames;  // (G) NumberOfGlobalVariables
799   float* GlobalVariableValue;  // (G) NumberOfGlobalVariables
800 
801   // The element and node arrays in the file were all scalar arrays.
802   // Those with similar names were combined into vectors in VTK.  Here
803   // are all the original names from the Exodus file, the names given
804   // the variables in the VTK ugrid, and a mapping from the VTK names
805   // to the Exodus names.
806 
807   int OriginalNumberOfElementVariables;   // (G)
808   char** OriginalElementVariableNames;    // (G) OriginalNumberOfElementVariables
809   int NumberOfElementVariables;           // (G)
810   int MaxNumberOfElementVariables;        // (G)
811   char** ElementVariableNames;            // (G) MaxNumberOfElementVariables
812   int* ElementVariableNumberOfComponents; // (G) MaxNumberOfElementVariables
813   int* MapToOriginalElementVariableNames; // (G) MaxNumberOfElementVariables
814 
815   int OriginalNumberOfNodeVariables;   // (G)
816   char** OriginalNodeVariableNames;    // (G) OriginalNumberOfNodeVariables
817   int NumberOfNodeVariables;           // (G)
818   int MaxNumberOfNodeVariables;        // (G)
819   char** NodeVariableNames;            // (G) NumberOfNodeVariables
820   int* NodeVariableNumberOfComponents; // (G) NumberOfNodeVariables
821   int* MapToOriginalNodeVariableNames; // (G) NumberOfNodeVariables
822 
823   int* ElementVariableTruthTable; // (G) NumBlocks*OrigNumberOfElementVariables
824   vtkTypeBool AllVariablesDefinedInAllBlocks;
825 
826 private:
827   vtkModelMetadata(const vtkModelMetadata&) = delete;
828   void operator=(const vtkModelMetadata&) = delete;
829 };
830 #endif
831