VTK  9.4.20250222
vtkLSDynaReader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
140#ifndef vtkLSDynaReader_h
141#define vtkLSDynaReader_h
142
143#include "vtkIOLSDynaModule.h" // For export macro
145#include <string> // for method signature
146
147VTK_ABI_NAMESPACE_BEGIN
148class LSDynaMetaData;
150class vtkPoints;
151class vtkDataArray;
154
155class VTKIOLSDYNA_EXPORT vtkLSDynaReader : public vtkMultiBlockDataSetAlgorithm
156{
157public:
159 void PrintSelf(ostream& os, vtkIndent indent) override;
161
166 void Dump(ostream& os);
167
172 void DebugDump();
173
177 virtual int CanReadFile(VTK_FILEPATH const char* fname);
178
180
184 virtual void SetDatabaseDirectory(VTK_FILEPATH const std::string&);
185 virtual void SetDatabaseDirectory(VTK_FILEPATH const char*);
189
191
197 virtual void SetFileName(VTK_FILEPATH const std::string&);
198 virtual void SetFileName(VTK_FILEPATH const char*);
201
207 char* GetTitle();
208
215
222
232
244
251
258
265
272
279
286
293
295
301 virtual void SetTimeStep(vtkIdType);
304 vtkGetVector2Macro(TimeStepRange, int);
305 vtkSetVector2Macro(TimeStepRange, int);
307
309
314 const char* GetPointArrayName(int);
315 virtual void SetPointArrayStatus(int arr, int status);
316 virtual void SetPointArrayStatus(const char* arrName, int status);
318 int GetPointArrayStatus(const char* arrName);
320 int GetNumberOfComponentsInPointArray(const char* arrName);
322
324
330 int GetNumberOfCellArrays(int cellType);
331 const char* GetCellArrayName(int cellType, int arr);
332 virtual void SetCellArrayStatus(int cellType, int arr, int status);
333 virtual void SetCellArrayStatus(int cellType, const char* arrName, int status);
334 int GetCellArrayStatus(int cellType, int arr);
335 int GetCellArrayStatus(int cellType, const char* arrName);
336 int GetNumberOfComponentsInCellArray(int cellType, int arr);
337 int GetNumberOfComponentsInCellArray(int cellType, const char* arrName);
339
341
346 const char* GetSolidArrayName(int);
347 virtual void SetSolidArrayStatus(int arr, int status);
348 virtual void SetSolidArrayStatus(const char* arrName, int status);
350 int GetSolidArrayStatus(const char* arrName);
352
354 int GetNumberOfComponentsInSolidArray(const char* arrName);
355
357
362 const char* GetThickShellArrayName(int);
363 virtual void SetThickShellArrayStatus(int arr, int status);
364 virtual void SetThickShellArrayStatus(const char* arrName, int status);
366 int GetThickShellArrayStatus(const char* arrName);
368
370 int GetNumberOfComponentsInThickShellArray(const char* arrName);
371
373
378 const char* GetShellArrayName(int);
379 virtual void SetShellArrayStatus(int arr, int status);
380 virtual void SetShellArrayStatus(const char* arrName, int status);
382 int GetShellArrayStatus(const char* arrName);
384
386 int GetNumberOfComponentsInShellArray(const char* arrName);
387
389
394 const char* GetRigidBodyArrayName(int);
395 virtual void SetRigidBodyArrayStatus(int arr, int status);
396 virtual void SetRigidBodyArrayStatus(const char* arrName, int status);
398 int GetRigidBodyArrayStatus(const char* arrName);
400
402 int GetNumberOfComponentsInRigidBodyArray(const char* arrName);
403
405
410 const char* GetRoadSurfaceArrayName(int);
411 virtual void SetRoadSurfaceArrayStatus(int arr, int status);
412 virtual void SetRoadSurfaceArrayStatus(const char* arrName, int status);
414 int GetRoadSurfaceArrayStatus(const char* arrName);
416
418 int GetNumberOfComponentsInRoadSurfaceArray(const char* arrName);
419
421
426 const char* GetBeamArrayName(int);
427 virtual void SetBeamArrayStatus(int arr, int status);
428 virtual void SetBeamArrayStatus(const char* arrName, int status);
429 int GetBeamArrayStatus(int arr);
430 int GetBeamArrayStatus(const char* arrName);
432
434 int GetNumberOfComponentsInBeamArray(const char* arrName);
435
437
442 const char* GetParticleArrayName(int);
443 virtual void SetParticleArrayStatus(int arr, int status);
444 virtual void SetParticleArrayStatus(const char* arrName, int status);
446 int GetParticleArrayStatus(const char* arrName);
448
450 int GetNumberOfComponentsInParticleArray(const char* arrName);
451
453
459 vtkGetMacro(DeformedMesh, vtkTypeBool);
460 vtkBooleanMacro(DeformedMesh, vtkTypeBool);
462
464
474 vtkSetMacro(RemoveDeletedCells, vtkTypeBool);
475 vtkGetMacro(RemoveDeletedCells, vtkTypeBool);
476 vtkBooleanMacro(RemoveDeletedCells, vtkTypeBool);
478
480
484 vtkSetMacro(DeletedCellsAsGhostArray, vtkTypeBool);
485 vtkGetMacro(DeletedCellsAsGhostArray, vtkTypeBool);
486 vtkBooleanMacro(DeletedCellsAsGhostArray, vtkTypeBool);
488
490
501 vtkSetStringMacro(InputDeck);
502 vtkGetStringMacro(InputDeck);
504
506
517 const char* GetPartArrayName(int);
518 virtual void SetPartArrayStatus(int arr, int status);
519 virtual void SetPartArrayStatus(const char* partName, int status);
520 int GetPartArrayStatus(int arr);
521 int GetPartArrayStatus(const char* partName);
523
524protected:
525 // holds all the parts and all the properties for each part
527
533
535
542
547 int TimeStepRange[2];
548
553
556
565 int ReadHeaderInformation(int currentAdaptLevel);
566
577
580
582
591 virtual int ReadTopology();
592 virtual int ReadNodes();
593 virtual int ReadPartSizes();
595 virtual int ReadUserIds();
596 virtual int ReadState(vtkIdType);
599 virtual int ReadDeletion();
603
607 virtual void ResetPartInfo();
608
613 virtual int ReadInputDeck();
614
621
627 virtual int ReadUserMaterialIds();
628
630
634 int ReadInputDeckXML(istream& deck);
635 int ReadInputDeckKeywords(istream& deck);
637
642 int WriteInputDeckSummary(const char* fname);
643
655 virtual void ReadDeletionArray(vtkUnsignedCharArray* arr, const int& pos, const int& size);
656
660 virtual void ReadCellProperties(const int& type, const int& numTuples);
661
663
665
666private:
667 // Helper templated methods to optimize reading. We cast the entire buffer
668 // to a given type instead of casting each element to improve performance
669 template <typename T>
670 void FillDeletionArray(T* buffer, vtkUnsignedCharArray* arr, const vtkIdType& start,
671 const vtkIdType& numCells, const int& deathPos, const int& cellSize);
672
673 template <int wordSize, typename T>
674 int FillTopology(T* buffer);
675
676 template <typename T, int blockType, vtkIdType numWordsPerCell, vtkIdType cellLength>
677 void ReadBlockCellSizes();
678
679 template <typename T>
680 int FillPartSizes();
681
682 vtkLSDynaReader(const vtkLSDynaReader&) = delete;
683 void operator=(const vtkLSDynaReader&) = delete;
684};
685
686inline void vtkLSDynaReader::SetPointArrayStatus(const char* arrName, int status)
687{
688 for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
689 {
690 if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
691 {
692 this->SetPointArrayStatus(a, status);
693 return;
694 }
695 }
696 vtkWarningMacro("Point array \"" << arrName << "\" does not exist");
697}
698
699inline int vtkLSDynaReader::GetPointArrayStatus(const char* arrName)
700{
701 for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
702 {
703 if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
704 {
705 return this->GetPointArrayStatus(a);
706 }
707 }
708 // vtkWarningMacro( "Point array \"" << arrName << "\" does not exist" );
709 return 0;
710}
711
713{
714 for (int a = 0; a < this->GetNumberOfPointArrays(); ++a)
715 {
716 if (strcmp(arrName, this->GetPointArrayName(a)) == 0)
717 {
718 return this->GetNumberOfComponentsInPointArray(a);
719 }
720 }
721 // vtkWarningMacro( "Point array \"" << arrName << "\" does not exist" );
722 return 0;
723}
724
725inline void vtkLSDynaReader::SetCellArrayStatus(int cellType, const char* arrName, int status)
726{
727 for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
728 {
729 if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
730 {
731 this->SetCellArrayStatus(cellType, a, status);
732 return;
733 }
734 }
735 vtkWarningMacro("Cell array \"" << arrName << "\" (type " << cellType << ") does not exist");
736}
737
738inline int vtkLSDynaReader::GetCellArrayStatus(int cellType, const char* arrName)
739{
740 for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
741 {
742 if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
743 {
744 return this->GetCellArrayStatus(cellType, a);
745 }
746 }
747 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
748 return 0;
749}
750
751inline int vtkLSDynaReader::GetNumberOfComponentsInCellArray(int cellType, const char* arrName)
752{
753 for (int a = 0; a < this->GetNumberOfCellArrays(cellType); ++a)
754 {
755 if (strcmp(arrName, this->GetCellArrayName(cellType, a)) == 0)
756 {
757 return this->GetNumberOfComponentsInCellArray(cellType, a);
758 }
759 }
760 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
761 return 0;
762}
763
764inline void vtkLSDynaReader::SetSolidArrayStatus(const char* arrName, int status)
765{
766 for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
767 {
768 if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
769 {
770 this->SetSolidArrayStatus(a, status);
771 return;
772 }
773 }
774 vtkWarningMacro("Solid array \"" << arrName << "\" does not exist");
775}
776
777inline int vtkLSDynaReader::GetSolidArrayStatus(const char* arrName)
778{
779 for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
780 {
781 if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
782 {
783 return this->GetSolidArrayStatus(a);
784 }
785 }
786 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
787 return 0;
788}
789
791{
792 for (int a = 0; a < this->GetNumberOfSolidArrays(); ++a)
793 {
794 if (strcmp(arrName, this->GetSolidArrayName(a)) == 0)
795 {
796 return this->GetNumberOfComponentsInSolidArray(a);
797 }
798 }
799 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
800 return 0;
801}
802
803inline void vtkLSDynaReader::SetThickShellArrayStatus(const char* arrName, int status)
804{
805 for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
806 {
807 if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
808 {
809 this->SetThickShellArrayStatus(a, status);
810 return;
811 }
812 }
813 vtkWarningMacro("Thick shell array \"" << arrName << "\" does not exist");
814}
815
816inline int vtkLSDynaReader::GetThickShellArrayStatus(const char* arrName)
817{
818 for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
819 {
820 if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
821 {
822 return this->GetThickShellArrayStatus(a);
823 }
824 }
825 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
826 return 0;
827}
828
830{
831 for (int a = 0; a < this->GetNumberOfThickShellArrays(); ++a)
832 {
833 if (strcmp(arrName, this->GetThickShellArrayName(a)) == 0)
834 {
836 }
837 }
838 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
839 return 0;
840}
841
842inline void vtkLSDynaReader::SetShellArrayStatus(const char* arrName, int status)
843{
844 for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
845 {
846 if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
847 {
848 this->SetShellArrayStatus(a, status);
849 return;
850 }
851 }
852 vtkWarningMacro("Shell array \"" << arrName << "\" does not exist");
853}
854
855inline int vtkLSDynaReader::GetShellArrayStatus(const char* arrName)
856{
857 for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
858 {
859 if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
860 {
861 return this->GetShellArrayStatus(a);
862 }
863 }
864 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
865 return 0;
866}
867
869{
870 for (int a = 0; a < this->GetNumberOfShellArrays(); ++a)
871 {
872 if (strcmp(arrName, this->GetShellArrayName(a)) == 0)
873 {
874 return this->GetNumberOfComponentsInShellArray(a);
875 }
876 }
877 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
878 return 0;
879}
880
881inline void vtkLSDynaReader::SetBeamArrayStatus(const char* arrName, int status)
882{
883 for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
884 {
885 if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
886 {
887 this->SetBeamArrayStatus(a, status);
888 return;
889 }
890 }
891 vtkWarningMacro("Beam array \"" << arrName << "\" does not exist");
892}
893
894inline int vtkLSDynaReader::GetBeamArrayStatus(const char* arrName)
895{
896 for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
897 {
898 if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
899 {
900 return this->GetBeamArrayStatus(a);
901 }
902 }
903 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
904 return 0;
905}
906
908{
909 for (int a = 0; a < this->GetNumberOfBeamArrays(); ++a)
910 {
911 if (strcmp(arrName, this->GetBeamArrayName(a)) == 0)
912 {
913 return this->GetNumberOfComponentsInBeamArray(a);
914 }
915 }
916 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
917 return 0;
918}
919
920inline void vtkLSDynaReader::SetParticleArrayStatus(const char* arrName, int status)
921{
922 for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
923 {
924 if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
925 {
926 this->SetParticleArrayStatus(a, status);
927 return;
928 }
929 }
930 vtkWarningMacro("Particle array \"" << arrName << "\" does not exist");
931}
932
933inline int vtkLSDynaReader::GetParticleArrayStatus(const char* arrName)
934{
935 for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
936 {
937 if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
938 {
939 return this->GetParticleArrayStatus(a);
940 }
941 }
942 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
943 return 0;
944}
945
947{
948 for (int a = 0; a < this->GetNumberOfParticleArrays(); ++a)
949 {
950 if (strcmp(arrName, this->GetParticleArrayName(a)) == 0)
951 {
953 }
954 }
955 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
956 return 0;
957}
958
959inline void vtkLSDynaReader::SetRigidBodyArrayStatus(const char* arrName, int status)
960{
961 for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
962 {
963 if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
964 {
965 this->SetRigidBodyArrayStatus(a, status);
966 return;
967 }
968 }
969 vtkWarningMacro("Rigid body array \"" << arrName << "\" does not exist");
970}
971
972inline int vtkLSDynaReader::GetRigidBodyArrayStatus(const char* arrName)
973{
974 for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
975 {
976 if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
977 {
978 return this->GetRigidBodyArrayStatus(a);
979 }
980 }
981 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
982 return 0;
983}
984
986{
987 for (int a = 0; a < this->GetNumberOfRigidBodyArrays(); ++a)
988 {
989 if (strcmp(arrName, this->GetRigidBodyArrayName(a)) == 0)
990 {
992 }
993 }
994 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
995 return 0;
996}
997
998inline void vtkLSDynaReader::SetRoadSurfaceArrayStatus(const char* arrName, int status)
999{
1000 for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1001 {
1002 if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1003 {
1004 this->SetRoadSurfaceArrayStatus(a, status);
1005 return;
1006 }
1007 }
1008 vtkWarningMacro("Road surface array \"" << arrName << "\" does not exist");
1009}
1010
1011inline int vtkLSDynaReader::GetRoadSurfaceArrayStatus(const char* arrName)
1012{
1013 for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1014 {
1015 if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1016 {
1017 return this->GetRoadSurfaceArrayStatus(a);
1018 }
1019 }
1020 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1021 return 0;
1022}
1023
1025{
1026 for (int a = 0; a < this->GetNumberOfRoadSurfaceArrays(); ++a)
1027 {
1028 if (strcmp(arrName, this->GetRoadSurfaceArrayName(a)) == 0)
1029 {
1031 }
1032 }
1033 // vtkWarningMacro( "Cell array \"" << arrName << "\" does not exist" );
1034 return 0;
1035}
1036
1037inline void vtkLSDynaReader::SetPartArrayStatus(const char* arrName, int status)
1038{
1039 for (int a = 0; a < this->GetNumberOfPartArrays(); ++a)
1040 {
1041 if (strcmp(arrName, this->GetPartArrayName(a)) == 0)
1042 {
1043 this->SetPartArrayStatus(a, status);
1044 return;
1045 }
1046 }
1047 vtkWarningMacro("Part \"" << arrName << "\" does not exist");
1048}
1049
1050inline int vtkLSDynaReader::GetPartArrayStatus(const char* partName)
1051{
1052 for (int a = 0; a < this->GetNumberOfPartArrays(); ++a)
1053 {
1054 if (strcmp(partName, this->GetPartArrayName(a)) == 0)
1055 {
1056 return this->GetPartArrayStatus(a);
1057 }
1058 }
1059 // vtkWarningMacro( "PartArray \"" << partName << "\" does not exist" );
1060 return 0;
1061}
1062
1063VTK_ABI_NAMESPACE_END
1064#endif // vtkLSDynaReader_h
abstract superclass for arrays of numeric data
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Read LS-Dyna databases (d3plot)
virtual void SetDatabaseDirectory(VTK_FILEPATH const char *)
Get/Set the directory containing the LS-Dyna database and determine whether it is valid.
const char * GetRigidBodyArrayName(int)
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.
int GetRigidBodyArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadSPHState(vtkIdType)
These functions read various parts of the database.
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...
virtual int ReadInputDeck()
Called from within ReadHeaderInformation() to read part names associated with material IDs.
virtual int ReadCellStateInfo(vtkIdType)
These functions read various parts of the database.
const char * GetSolidArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
void Dump(ostream &os)
Print out more complete information about the dataset (and less complete information about the VTK hi...
int GetRoadSurfaceArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadPartTitlesFromRootFile()
Called from within ReadHeaderInformation to read part names from the end of the first d3plot file.
vtkIdType GetNumberOfCells()
Retrieve the number of cells of a given type in the database.
int GetNumberOfComponentsInParticleArray(int a)
int GetParticleArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadUserIds()
These functions read various parts of the database.
virtual void SetThickShellArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
vtkIdType GetNumberOfShellCells()
Retrieve the number of cells of a given type in the database.
int ReadInputDeckKeywords(istream &deck)
ReadInputDeck determines the type of file (keyword or XML summary) and calls one of these two routine...
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 SetShellArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
int GetNumberOfComponentsInBeamArray(int a)
int ReadHeaderInformation(int currentAdaptLevel)
This function populates the reader's private dictionary with information about the database.
virtual void SetDatabaseDirectory(VTK_FILEPATH const std::string &)
Get/Set the directory containing the LS-Dyna database and determine whether it is valid.
const char * GetRoadSurfaceArrayName(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 int CanReadFile(VTK_FILEPATH const char *fname)
Determine if the file can be read with this reader.
virtual void SetSolidArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
void ResetPartsCache()
virtual int ReadUserMaterialIds()
Called from within ReadHeaderInformation() to read arbitrary material IDs (if present) or manufacture...
int GetThickShellArrayStatus(int arr)
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.
char * GetTitle()
The title of the database is a 40 or 80 character text description stored at the front of a d3plot fi...
int GetBeamArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadNodeStateInfo(vtkIdType)
These functions read various parts of the database.
vtkLSDynaPartCollection * Parts
void DebugDump()
A routine to call Dump() from within a lame debugger that won't properly pass a C++ iostream object l...
int GetNumberOfRigidBodyArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadState(vtkIdType)
These functions read various parts of the database.
int GetNumberOfPartArrays()
These methods allow you to load only selected parts of the input.
const char * GetThickShellArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
const char * GetParticleArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
vtkIdType GetNumberOfRigidBodyCells()
Retrieve the number of cells of a given type in the database.
int GetNumberOfComponentsInPointArray(int arr)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh.
vtkTypeBool DeformedMesh
Should deflected coordinates be used, or should the mesh remain undeflected? By default,...
const char * GetPointArrayName(int)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh.
virtual void ReadDeletionArray(vtkUnsignedCharArray *arr, const int &pos, const int &size)
Read an array of deletion data.
int WriteInputDeckSummary(const char *fname)
ReadInputDeckKeywords calls this function if it was successful in reading part names for materials.
int GetDimensionality()
Retrieve the dimension of points in the database.
const char * GetShellArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
vtkIdType GetNumberOfThickShellCells()
Retrieve the number of cells of a given type in the database.
virtual void SetParticleArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
static vtkLSDynaReader * New()
int GetNumberOfShellArrays()
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.
~vtkLSDynaReader() override
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 GetNumberOfComponentsInRigidBodyArray(int a)
virtual void SetFileName(VTK_FILEPATH const char *)
Get/Set the filename.
const char * GetPartArrayName(int)
These methods allow you to load only selected parts of the input.
int GetNumberOfThickShellArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
double GetTimeValue(vtkIdType)
Retrieve information about the time extents of the LS-Dyna database.
int GetNumberOfCellArrays(int cellType)
Routines that allow the status of a cell variable to be adjusted or queried independent of the output...
int GetNumberOfRoadSurfaceArrays()
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual int ReadPartSizes()
These functions read various parts of the database.
virtual void SetBeamArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
int ReadInputDeckXML(istream &deck)
ReadInputDeck determines the type of file (keyword or XML summary) and calls one of these two routine...
vtkIdType GetTimeStep()
Retrieve information about the time extents of the LS-Dyna database.
int IsDatabaseValid()
Get/Set the directory containing the LS-Dyna database and determine whether it is valid.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIdType GetNumberOfParticleCells()
Retrieve the number of cells of a given type in the database.
VTK_FILEPATH std::string GetDatabaseDirectory()
Get/Set the directory containing the LS-Dyna database and determine whether it is valid.
int GetNumberOfBeamArrays()
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.
vtkTypeBool DeletedCellsAsGhostArray
Should cells marked as deleted be removed from the mesh? By default, this is true.
virtual int ReadConnectivityAndMaterial()
These functions read various parts of the database.
virtual void ResetPartInfo()
Resets the Part information to the default state.
virtual void SetFileName(VTK_FILEPATH const std::string &)
Get/Set the filename.
virtual int ReadNodes()
These functions read various parts of the database.
vtkIdType GetNumberOfNodes()
Retrieve the number of points in the database.
virtual void SetRoadSurfaceArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual int ReadTopology()
These functions read various parts of the database.
VTK_FILEPATH std::string GetFileName()
Get/Set the filename.
int GetNumberOfComponentsInRoadSurfaceArray(int a)
int GetNumberOfComponentsInSolidArray(int a)
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int GetNumberOfPointArrays()
These methods allow you to load only selected subsets of the nodal variables defined over the mesh.
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...
int GetNumberOfComponentsInThickShellArray(int a)
virtual int ComputeDeflectionAndUpdateGeometry(vtkUnstructuredGrid *grid)
These functions read various parts of the database.
vtkIdType GetNumberOfSolidCells()
Retrieve the number of cells of a given type in the database.
LSDynaMetaData * P
virtual void ReadCellProperties(const int &type, const int &numTuples)
Read all the cell properties of a given part type.
int ScanDatabaseTimeSteps()
This function scans the list of files in the database and bookmarks the start of each time step's sta...
vtkIdType GetNumberOfTimeSteps()
Retrieve information about the time extents of the LS-Dyna database.
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 GetPointArrayStatus(int arr)
These methods allow you to load only selected subsets of the nodal variables defined over the mesh.
vtkIdType GetNumberOfRoadSurfaceCells()
Retrieve the number of cells of a given type in the database.
vtkIdType GetNumberOfBeamCells()
Retrieve the number of cells of a given type in the database.
int GetNumberOfComponentsInShellArray(int a)
virtual void SetTimeStep(vtkIdType)
Retrieve information about the time extents of the LS-Dyna database.
vtkTypeBool RemoveDeletedCells
Should cells marked as deleted be removed from the mesh? By default, this is true.
int GetSolidArrayStatus(int arr)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
virtual void SetRigidBodyArrayStatus(int arr, int status)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
vtkIdType GetNumberOfContinuumCells()
Retrieve the number of cells of a given type in the database.
void SetDeformedMesh(vtkTypeBool)
Should deflected coordinates be used, or should the mesh remain undeflected? By default,...
virtual int ReadDeletion()
These functions read various parts of the database.
int GetPartArrayStatus(int arr)
These methods allow you to load only selected parts of the input.
const char * GetBeamArrayName(int)
These methods allow you to load only selected subsets of the cell variables defined over the mesh.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
represent and manipulate 3D points
Definition vtkPoints.h:139
dynamic, self-adjusting array of unsigned char
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
#define VTK_FILEPATH