VTK  9.0.1
vtkLSDynaReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLSDynaReader.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 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
152 #ifndef vtkLSDynaReader_h
153 #define vtkLSDynaReader_h
154 
155 #include "vtkIOLSDynaModule.h" // For export macro
157 #include <string> // for method signature
158 
159 class LSDynaMetaData;
161 class vtkPoints;
162 class vtkDataArray;
164 class vtkUnstructuredGrid;
165 
166 class VTKIOLSDYNA_EXPORT vtkLSDynaReader : public vtkMultiBlockDataSetAlgorithm
167 {
168 public:
170  void PrintSelf(ostream& os, vtkIndent indent) override;
171  static vtkLSDynaReader* New();
172 
177  void Dump(ostream& os);
178 
183  void DebugDump();
184 
188  virtual int CanReadFile(const char* fname);
189 
191 
195  virtual void SetDatabaseDirectory(const std::string&);
196  virtual void SetDatabaseDirectory(const char*);
197 #ifdef VTK_LEGACY_REMOVE
198  std::string GetDatabaseDirectory();
199 #else
200  VTK_LEGACY(const char* GetDatabaseDirectory());
201 #endif
202  int IsDatabaseValid();
204 
206 
212  virtual void SetFileName(const std::string&);
213  virtual void SetFileName(const char*);
214 #ifdef VTK_LEGACY_REMOVE
216 #else
217  VTK_LEGACY(const char* GetFileName());
218 #endif
219 
220 
226  char* GetTitle();
227 
233  int GetDimensionality();
234 
240  vtkIdType GetNumberOfNodes();
241 
250  vtkIdType GetNumberOfCells();
251 
262  vtkIdType GetNumberOfContinuumCells();
263 
269  vtkIdType GetNumberOfSolidCells();
270 
276  vtkIdType GetNumberOfThickShellCells();
277 
283  vtkIdType GetNumberOfShellCells();
284 
290  vtkIdType GetNumberOfRigidBodyCells();
291 
297  vtkIdType GetNumberOfRoadSurfaceCells();
298 
304  vtkIdType GetNumberOfBeamCells();
305 
311  vtkIdType GetNumberOfParticleCells();
312 
314 
319  vtkIdType GetNumberOfTimeSteps();
320  virtual void SetTimeStep(vtkIdType);
321  vtkIdType GetTimeStep();
322  double GetTimeValue(vtkIdType);
323  vtkGetVector2Macro(TimeStepRange, int);
324  vtkSetVector2Macro(TimeStepRange, int);
326 
328 
332  int GetNumberOfPointArrays();
333  const char* GetPointArrayName(int);
334  virtual void SetPointArrayStatus(int arr, int status);
335  virtual void SetPointArrayStatus(const char* arrName, int status);
336  int GetPointArrayStatus(int arr);
337  int GetPointArrayStatus(const char* arrName);
338  int GetNumberOfComponentsInPointArray(int arr);
339  int GetNumberOfComponentsInPointArray(const char* arrName);
341 
343 
349  int GetNumberOfCellArrays(int cellType);
350  const char* GetCellArrayName(int cellType, int arr);
351  virtual void SetCellArrayStatus(int cellType, int arr, int status);
352  virtual void SetCellArrayStatus(int cellType, const char* arrName, int status);
353  int GetCellArrayStatus(int cellType, int arr);
354  int GetCellArrayStatus(int cellType, const char* arrName);
355  int GetNumberOfComponentsInCellArray(int cellType, int arr);
356  int GetNumberOfComponentsInCellArray(int cellType, const char* arrName);
358 
360 
364  int GetNumberOfSolidArrays();
365  const char* GetSolidArrayName(int);
366  virtual void SetSolidArrayStatus(int arr, int status);
367  virtual void SetSolidArrayStatus(const char* arrName, int status);
368  int GetSolidArrayStatus(int arr);
369  int GetSolidArrayStatus(const char* arrName);
371 
372  int GetNumberOfComponentsInSolidArray(int a);
373  int GetNumberOfComponentsInSolidArray(const char* arrName);
374 
376 
380  int GetNumberOfThickShellArrays();
381  const char* GetThickShellArrayName(int);
382  virtual void SetThickShellArrayStatus(int arr, int status);
383  virtual void SetThickShellArrayStatus(const char* arrName, int status);
384  int GetThickShellArrayStatus(int arr);
385  int GetThickShellArrayStatus(const char* arrName);
387 
388  int GetNumberOfComponentsInThickShellArray(int a);
389  int GetNumberOfComponentsInThickShellArray(const char* arrName);
390 
392 
396  int GetNumberOfShellArrays();
397  const char* GetShellArrayName(int);
398  virtual void SetShellArrayStatus(int arr, int status);
399  virtual void SetShellArrayStatus(const char* arrName, int status);
400  int GetShellArrayStatus(int arr);
401  int GetShellArrayStatus(const char* arrName);
403 
404  int GetNumberOfComponentsInShellArray(int a);
405  int GetNumberOfComponentsInShellArray(const char* arrName);
406 
408 
412  int GetNumberOfRigidBodyArrays();
413  const char* GetRigidBodyArrayName(int);
414  virtual void SetRigidBodyArrayStatus(int arr, int status);
415  virtual void SetRigidBodyArrayStatus(const char* arrName, int status);
416  int GetRigidBodyArrayStatus(int arr);
417  int GetRigidBodyArrayStatus(const char* arrName);
419 
420  int GetNumberOfComponentsInRigidBodyArray(int a);
421  int GetNumberOfComponentsInRigidBodyArray(const char* arrName);
422 
424 
428  int GetNumberOfRoadSurfaceArrays();
429  const char* GetRoadSurfaceArrayName(int);
430  virtual void SetRoadSurfaceArrayStatus(int arr, int status);
431  virtual void SetRoadSurfaceArrayStatus(const char* arrName, int status);
432  int GetRoadSurfaceArrayStatus(int arr);
433  int GetRoadSurfaceArrayStatus(const char* arrName);
435 
436  int GetNumberOfComponentsInRoadSurfaceArray(int a);
437  int GetNumberOfComponentsInRoadSurfaceArray(const char* arrName);
438 
440 
444  int GetNumberOfBeamArrays();
445  const char* GetBeamArrayName(int);
446  virtual void SetBeamArrayStatus(int arr, int status);
447  virtual void SetBeamArrayStatus(const char* arrName, int status);
448  int GetBeamArrayStatus(int arr);
449  int GetBeamArrayStatus(const char* arrName);
451 
452  int GetNumberOfComponentsInBeamArray(int a);
453  int GetNumberOfComponentsInBeamArray(const char* arrName);
454 
456 
460  int GetNumberOfParticleArrays();
461  const char* GetParticleArrayName(int);
462  virtual void SetParticleArrayStatus(int arr, int status);
463  virtual void SetParticleArrayStatus(const char* arrName, int status);
464  int GetParticleArrayStatus(int arr);
465  int GetParticleArrayStatus(const char* arrName);
467 
468  int GetNumberOfComponentsInParticleArray(int a);
469  int GetNumberOfComponentsInParticleArray(const char* arrName);
470 
472 
477  void SetDeformedMesh(vtkTypeBool);
478  vtkGetMacro(DeformedMesh, vtkTypeBool);
479  vtkBooleanMacro(DeformedMesh, vtkTypeBool);
481 
483 
493  vtkSetMacro(RemoveDeletedCells, vtkTypeBool);
494  vtkGetMacro(RemoveDeletedCells, vtkTypeBool);
495  vtkBooleanMacro(RemoveDeletedCells, vtkTypeBool);
497 
499 
503  vtkSetMacro(DeletedCellsAsGhostArray, vtkTypeBool);
504  vtkGetMacro(DeletedCellsAsGhostArray, vtkTypeBool);
505  vtkBooleanMacro(DeletedCellsAsGhostArray, vtkTypeBool);
507 
509 
520  vtkSetStringMacro(InputDeck);
521  vtkGetStringMacro(InputDeck);
523 
525 
535  int GetNumberOfPartArrays();
536  const char* GetPartArrayName(int);
537  virtual void SetPartArrayStatus(int arr, int status);
538  virtual void SetPartArrayStatus(const char* partName, int status);
539  int GetPartArrayStatus(int arr);
540  int GetPartArrayStatus(const char* partName);
542 
543 protected:
544  // holds all the parts and all the properties for each part
546 
552 
554 
561 
566  int TimeStepRange[2];
567 
571  char* InputDeck;
572 
573  vtkLSDynaReader();
574  ~vtkLSDynaReader() override;
575 
584  int ReadHeaderInformation(int currentAdaptLevel);
585 
595  int ScanDatabaseTimeSteps();
596 
599 
601 
610  virtual int ReadTopology();
611  virtual int ReadNodes();
612  virtual int ReadPartSizes();
613  virtual int ReadConnectivityAndMaterial();
614  virtual int ReadUserIds();
615  virtual int ReadState(vtkIdType);
616  virtual int ReadNodeStateInfo(vtkIdType);
617  virtual int ReadCellStateInfo(vtkIdType);
618  virtual int ReadDeletion();
619  virtual int ReadSPHState(vtkIdType);
620  virtual int ComputeDeflectionAndUpdateGeometry(vtkUnstructuredGrid* grid);
622 
626  virtual void ResetPartInfo();
627 
632  virtual int ReadInputDeck();
633 
639  virtual int ReadPartTitlesFromRootFile();
640 
646  virtual int ReadUserMaterialIds();
647 
649 
653  int ReadInputDeckXML(istream& deck);
654  int ReadInputDeckKeywords(istream& deck);
656 
661  int WriteInputDeckSummary(const char* fname);
662 
674  virtual void ReadDeletionArray(vtkUnsignedCharArray* arr, const int& pos, const int& size);
675 
679  virtual void ReadCellProperties(const int& type, const int& numTuples);
680 
682 
683  void ResetPartsCache();
684 
685 private:
686  // Helper templated methods to optimize reading. We cast the entire buffer
687  // to a given type instead of casting each element to improve performance
688  template <typename T>
689  void FillDeletionArray(T* buffer, vtkUnsignedCharArray* arr, const vtkIdType& start,
690  const vtkIdType& numCells, const int& deathPos, const int& cellSize);
691 
692  template <int wordSize, typename T>
693  int FillTopology(T* buffer);
694 
695  template <typename T, int blockType, vtkIdType numWordsPerCell, vtkIdType cellLength>
696  void ReadBlockCellSizes();
697 
698  template <typename T>
699  int FillPartSizes();
700 
701  vtkLSDynaReader(const vtkLSDynaReader&) = delete;
702  void operator=(const vtkLSDynaReader&) = delete;
703 };
704 
705 inline void vtkLSDynaReader::SetPointArrayStatus(const char* arrName, int status)
706 {
707  for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
708  {
709  if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
710  {
711  this->SetPointArrayStatus(a, status);
712  return;
713  }
714  }
715  vtkWarningMacro("Point array \"" << arrName << "\" does not exist");
716 }
717 
718 inline int vtkLSDynaReader::GetPointArrayStatus(const char* arrName)
719 {
720  for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
721  {
722  if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
723  {
724  return this->GetPointArrayStatus(a);
725  }
726  }
727  // vtkWarningMacro( "Point array \"" << arrName << "\" does not exist" );
728  return 0;
729 }
730 
732 {
733  for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
734  {
735  if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
736  {
737  return this->GetNumberOfComponentsInPointArray(a);
738  }
739  }
740  // vtkWarningMacro( "Point array \"" << arrName << "\" does not exist" );
741  return 0;
742 }
743 
744 inline void vtkLSDynaReader::SetCellArrayStatus(int cellType, const char* arrName, int status)
745 {
746  for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
747  {
748  if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
749  {
750  this->SetCellArrayStatus(cellType, a, status);
751  return;
752  }
753  }
754  vtkWarningMacro("Cell array \"" << arrName << "\" (type " << cellType << ") does not exist");
755 }
756 
757 inline int vtkLSDynaReader::GetCellArrayStatus(int cellType, const char* arrName)
758 {
759  for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
760  {
761  if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
762  {
763  return this->GetCellArrayStatus(cellType, a);
764  }
765  }
766  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
767  return 0;
768 }
769 
770 inline int vtkLSDynaReader::GetNumberOfComponentsInCellArray(int cellType, const char* arrName)
771 {
772  for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
773  {
774  if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
775  {
776  return this->GetNumberOfComponentsInCellArray(cellType, a);
777  }
778  }
779  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
780  return 0;
781 }
782 
783 inline void vtkLSDynaReader::SetSolidArrayStatus(const char* arrName, int status)
784 {
785  for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
786  {
787  if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
788  {
789  this->SetSolidArrayStatus(a, status);
790  return;
791  }
792  }
793  vtkWarningMacro("Solid array \"" << arrName << "\" does not exist");
794 }
795 
796 inline int vtkLSDynaReader::GetSolidArrayStatus(const char* arrName)
797 {
798  for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
799  {
800  if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
801  {
802  return this->GetSolidArrayStatus(a);
803  }
804  }
805  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
806  return 0;
807 }
808 
810 {
811  for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
812  {
813  if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
814  {
815  return this->GetNumberOfComponentsInSolidArray(a);
816  }
817  }
818  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
819  return 0;
820 }
821 
822 inline void vtkLSDynaReader::SetThickShellArrayStatus(const char* arrName, int status)
823 {
824  for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
825  {
826  if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
827  {
828  this->SetThickShellArrayStatus(a, status);
829  return;
830  }
831  }
832  vtkWarningMacro("Thick shell array \"" << arrName << "\" does not exist");
833 }
834 
835 inline int vtkLSDynaReader::GetThickShellArrayStatus(const char* arrName)
836 {
837  for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
838  {
839  if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
840  {
841  return this->GetThickShellArrayStatus(a);
842  }
843  }
844  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
845  return 0;
846 }
847 
849 {
850  for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
851  {
852  if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
853  {
855  }
856  }
857  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
858  return 0;
859 }
860 
861 inline void vtkLSDynaReader::SetShellArrayStatus(const char* arrName, int status)
862 {
863  for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
864  {
865  if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
866  {
867  this->SetShellArrayStatus(a, status);
868  return;
869  }
870  }
871  vtkWarningMacro("Shell array \"" << arrName << "\" does not exist");
872 }
873 
874 inline int vtkLSDynaReader::GetShellArrayStatus(const char* arrName)
875 {
876  for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
877  {
878  if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
879  {
880  return this->GetShellArrayStatus(a);
881  }
882  }
883  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
884  return 0;
885 }
886 
888 {
889  for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
890  {
891  if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
892  {
893  return this->GetNumberOfComponentsInShellArray(a);
894  }
895  }
896  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
897  return 0;
898 }
899 
900 inline void vtkLSDynaReader::SetBeamArrayStatus(const char* arrName, int status)
901 {
902  for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
903  {
904  if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
905  {
906  this->SetBeamArrayStatus(a, status);
907  return;
908  }
909  }
910  vtkWarningMacro("Beam array \"" << arrName << "\" does not exist");
911 }
912 
913 inline int vtkLSDynaReader::GetBeamArrayStatus(const char* arrName)
914 {
915  for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
916  {
917  if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
918  {
919  return this->GetBeamArrayStatus(a);
920  }
921  }
922  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
923  return 0;
924 }
925 
927 {
928  for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
929  {
930  if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
931  {
932  return this->GetNumberOfComponentsInBeamArray(a);
933  }
934  }
935  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
936  return 0;
937 }
938 
939 inline void vtkLSDynaReader::SetParticleArrayStatus(const char* arrName, int status)
940 {
941  for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
942  {
943  if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
944  {
945  this->SetParticleArrayStatus(a, status);
946  return;
947  }
948  }
949  vtkWarningMacro("Particle array \"" << arrName << "\" does not exist");
950 }
951 
952 inline int vtkLSDynaReader::GetParticleArrayStatus(const char* arrName)
953 {
954  for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
955  {
956  if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
957  {
958  return this->GetParticleArrayStatus(a);
959  }
960  }
961  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
962  return 0;
963 }
964 
966 {
967  for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
968  {
969  if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
970  {
971  return this->GetNumberOfComponentsInParticleArray(a);
972  }
973  }
974  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
975  return 0;
976 }
977 
978 inline void vtkLSDynaReader::SetRigidBodyArrayStatus(const char* arrName, int status)
979 {
980  for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
981  {
982  if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
983  {
984  this->SetRigidBodyArrayStatus(a, status);
985  return;
986  }
987  }
988  vtkWarningMacro("Rigid body array \"" << arrName << "\" does not exist");
989 }
990 
991 inline int vtkLSDynaReader::GetRigidBodyArrayStatus(const char* arrName)
992 {
993  for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
994  {
995  if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
996  {
997  return this->GetRigidBodyArrayStatus(a);
998  }
999  }
1000  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1001  return 0;
1002 }
1003 
1005 {
1006  for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
1007  {
1008  if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
1009  {
1010  return this->GetNumberOfComponentsInRigidBodyArray(a);
1011  }
1012  }
1013  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1014  return 0;
1015 }
1016 
1017 inline void vtkLSDynaReader::SetRoadSurfaceArrayStatus(const char* arrName, int status)
1018 {
1019  for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1020  {
1021  if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1022  {
1023  this->SetRoadSurfaceArrayStatus(a, status);
1024  return;
1025  }
1026  }
1027  vtkWarningMacro("Road surface array \"" << arrName << "\" does not exist");
1028 }
1029 
1030 inline int vtkLSDynaReader::GetRoadSurfaceArrayStatus(const char* arrName)
1031 {
1032  for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1033  {
1034  if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1035  {
1036  return this->GetRoadSurfaceArrayStatus(a);
1037  }
1038  }
1039  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1040  return 0;
1041 }
1042 
1044 {
1045  for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1046  {
1047  if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1048  {
1050  }
1051  }
1052  // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1053  return 0;
1054 }
1055 
1056 inline void vtkLSDynaReader::SetPartArrayStatus(const char* arrName, int status)
1057 {
1058  for (int a = 0; a < this->GetNumberOfPartArrays(); ++a)
1059  {
1060  if (strcmp(arrName, this->GetPartArrayName(a)) == 0)
1061  {
1062  this->SetPartArrayStatus(a, status);
1063  return;
1064  }
1065  }
1066  vtkWarningMacro("Part \"" << arrName << "\" does not exist");
1067 }
1068 
1069 inline int vtkLSDynaReader::GetPartArrayStatus(const char* partName)
1070 {
1071  for (int a = 0; a < this->GetNumberOfPartArrays(); ++a)
1072  {
1073  if (strcmp(partName, this->GetPartArrayName(a)) == 0)
1074  {
1075  return this->GetPartArrayStatus(a);
1076  }
1077  }
1078  // vtkWarningMacro( "PartArray \"" << partName << "\" does not exist" );
1079  return 0;
1080 }
1081 
1082 #endif // vtkLSDynaReader_h
virtual void SetThickShellArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetShellArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
vtkLSDynaPartCollection * Parts
const char * GetPartArrayName(int)
These methods allow you to load only selected parts of the input.
int GetNumberOfCellArrays(int cellType)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
int GetNumberOfPointArrays()
These methods allow you to load only selected subsets of the nodal variables defined over the mesh...
std::string GetFileName(const std::string &fileName) noexcept
Set the appropriate file name based on recognized user input.
int GetThickShellArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
const char * GetPointArrayName(int)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh...
Store vtkAlgorithm input/output information.
int GetNumberOfComponentsInBeamArray(int a)
int GetNumberOfBeamArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfSolidArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetRoadSurfaceArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
Read LS-Dyna databases (d3plot)
vtkTypeBool DeformedMesh
Should deflected coordinates be used, or should the mesh remain undeflected? By default, this is true.
const char * GetParticleArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
const char * GetThickShellArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfComponentsInCellArray(int cellType, int arr)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
virtual void SetRigidBodyArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfComponentsInSolidArray(int a)
int GetNumberOfThickShellArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
const char * GetSolidArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int vtkIdType
Definition: vtkType.h:338
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
virtual void SetSolidArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
vtkTypeBool DeletedCellsAsGhostArray
Should cells marked as deleted be removed from the mesh? By default, this is true.
virtual void SetRoadSurfaceArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
static vtkMultiBlockDataSetAlgorithm * New()
int vtkTypeBool
Definition: vtkABI.h:69
int GetPointArrayStatus(int arr)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh...
int GetRigidBodyArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfParticleArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
virtual void SetShellArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetCellArrayStatus(int cellType, int arr)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
int GetNumberOfPartArrays()
These methods allow you to load only selected parts of the input.
int GetNumberOfRoadSurfaceArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
const char * GetRigidBodyArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
a simple class to control print indentation
Definition: vtkIndent.h:33
int GetNumberOfRigidBodyArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetSolidArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
int GetNumberOfShellArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkTypeBool RemoveDeletedCells
Should cells marked as deleted be removed from the mesh? By default, this is true.
dynamic, self-adjusting array of unsigned char
int GetNumberOfComponentsInShellArray(int a)
int GetNumberOfComponentsInRigidBodyArray(int a)
const char * GetCellArrayName(int cellType, int arr)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
int GetPartArrayStatus(int arr)
These methods allow you to load only selected parts of the input.
int GetParticleArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfComponentsInThickShellArray(int a)
int GetBeamArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
virtual void SetBeamArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
Store zero or more vtkInformation instances.
virtual void SetCellArrayStatus(int cellType, int arr, int status)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
const char * GetShellArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
char * InputDeck
The name of a file containing part names and IDs.
virtual void SetParticleArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
virtual void SetPartArrayStatus(int arr, int status)
These methods allow you to load only selected parts of the input.
virtual void SetPointArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh...
int GetNumberOfComponentsInRoadSurfaceArray(int a)
const char * GetBeamArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfComponentsInParticleArray(int a)
represent and manipulate 3D points
Definition: vtkPoints.h:33
const char * GetRoadSurfaceArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh...
int GetNumberOfComponentsInPointArray(int arr)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh...
virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
LSDynaMetaData * P