vtkbone
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkboneFiniteElementModel Class Reference

data model for finite element meshes More...

#include <vtkboneFiniteElementModel.h>

Inheritance diagram for vtkboneFiniteElementModel:
Inheritance graph
[legend]
Collaboration diagram for vtkboneFiniteElementModel:
Collaboration graph
[legend]

Public Types

enum  ElementType_t {
  N88_UNKNOWN, N88_TETRAHEDRON, N88_HEXAHEDRON, N88_MIXED,
  NUMBER_OF_ElementType
}
 
- Public Types inherited from vtkUnstructuredGrid
typedef vtkUnstructuredGridBase Superclass
 
- Public Types inherited from vtkPointSet
typedef vtkDataSet Superclass
 
- Public Types inherited from vtkDataSet
enum  FieldDataType
 
typedef vtkDataObject Superclass
 
- Public Types inherited from vtkDataObject
enum  FieldAssociations
 
enum  AttributeTypes
 
enum  FieldOperations
 
enum  FieldAssociations
 
enum  AttributeTypes
 
enum  FieldOperations
 
typedef vtkObject Superclass
 

Public Member Functions

virtual const char * GetClassName ()
 
virtual int IsA (const char *type)
 
void PrintSelf (ostream &os, vtkIndent indent)
 
virtual void AddNodeSet (vtkIdTypeArray *nodeSet)
 
virtual void AddElementSet (vtkIdTypeArray *elementSet)
 
virtual vtkIdTypeArrayGetNodeSet (const char *nodeSetName)
 
virtual vtkIdTypeArrayGetElementSet (const char *elementSetName)
 
virtual int RemoveNodeSet (const char *nodeSetName)
 
virtual int RemoveElementSet (const char *elementSetName)
 
virtual int GetElementType ()
 
virtual int ConvergenceSetFromConstraint (vtkboneConstraint *constraint)
 
void AppendHistory (const char *s)
 
void AppendLog (const char *s)
 
virtual int GetAssociatedElementsFromNodeSet (const char *nodeSetName, vtkIdTypeArray *ids)
 
virtual vtkIdTypeArrayGetAssociatedElementsFromNodeSet (const char *nodeSetName)
 
virtual int DataSetFromNodeSet (const vtkIdTypeArray *nodeSet, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromNodeSet (const vtkIdTypeArray *nodeSet)
 
virtual int DataSetFromNodeSet (const char *nodeSetName, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromNodeSet (const char *nodeSetName)
 
virtual int DataSetFromElementSet (const vtkIdTypeArray *elementSet, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromElementSet (const vtkIdTypeArray *elementSet)
 
virtual int DataSetFromElementSet (const char *elementSetName, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromElementSet (const char *elementSetName)
 
virtual void SetConstraints (vtkboneConstraintCollection *)
 
virtual vtkboneConstraintCollectionGetConstraints ()
 
virtual void SetMaterialTable (vtkboneMaterialTable *)
 
virtual vtkboneMaterialTableGetMaterialTable ()
 
virtual void SetConvergenceSet (vtkboneConstraint *)
 
virtual vtkboneConstraintGetConvergenceSet ()
 
virtual void GetAllCellPoints (vtkIdTypeArray *allCellPoints)
 
virtual vtkIdTypeArrayGetAllCellPoints ()
 
virtual int ApplyBoundaryCondition (vtkIdTypeArray *nodeIds, vtkDataArray *senses, vtkDataArray *displacements, const char *arg_constraintName)
 
virtual int ApplyBoundaryCondition (vtkIdTypeArray *nodeIds, int sense, double displacement, const char *arg_constraintName)
 
virtual int ApplyBoundaryCondition (vtkIdType nodeId, int sense, double displacement, const char *arg_constraintName)
 
virtual int ApplyBoundaryCondition (const char *nodeSetName, int sense, double displacement, const char *arg_constraintName)
 
virtual int FixNodes (vtkIdTypeArray *ids, const char *arg_constraintName)
 
virtual int FixNodes (vtkIdType id, const char *arg_constraintName)
 
virtual int FixNodes (const char *selectionName, const char *arg_constraintName)
 
virtual int ApplyLoad (vtkIdTypeArray *elementIds, vtkDataArray *distributions, vtkDataArray *senses, vtkDataArray *forces, const char *arg_constraintName)
 
virtual int ApplyLoad (vtkIdTypeArray *elementIds, int distribution, vtkDataArray *senses, vtkDataArray *forces, const char *arg_constraintName)
 
virtual int ApplyLoad (vtkIdTypeArray *elementIds, int distribution, int sense, double total_force, const char *arg_constraintName)
 
virtual int ApplyLoad (vtkIdType elementId, int distribution, int sense, double force, const char *arg_constraintName)
 
virtual int ApplyLoad (const char *elementSetName, int distribution, int sense, double total_force, const char *arg_constraintName)
 
virtual int DataSetFromConstraint (vtkboneConstraint *arg_constraint, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromConstraint (vtkboneConstraint *arg_constraint)
 
virtual int DataSetFromConstraint (const char *arg_constraintName, vtkUnstructuredGrid *data)
 
virtual vtkUnstructuredGridDataSetFromConstraint (const char *arg_constraintName)
 
virtual void SetName (const char *)
 
virtual char * GetName ()
 
virtual void SetNodeSets (vtkDataArrayCollection *)
 
virtual vtkDataArrayCollectionGetNodeSets ()
 
virtual void SetElementSets (vtkDataArrayCollection *)
 
virtual vtkDataArrayCollectionGetElementSets ()
 
virtual void SetGaussPointData (vtkDataArrayCollection *)
 
virtual vtkDataArrayCollectionGetGaussPointData ()
 
virtual void SetHistory (const char *)
 
virtual char * GetHistory ()
 
virtual void SetLog (const char *)
 
virtual char * GetLog ()
 
virtual void ShallowCopy (vtkDataObject *src)
 
virtual void DeepCopy (vtkDataObject *src)
 
- Public Member Functions inherited from vtkUnstructuredGrid
vtkUnstructuredGridNewInstance () const
 
int GetDataObjectType ()
 
virtual void Allocate (vtkIdType numCells=1000, int extSize=1000)
 
vtkIdType InsertNextCell (int type, vtkIdType npts, vtkIdType *ptIds)
 
vtkIdType InsertNextCell (int type, vtkIdList *ptIds)
 
vtkIdType InsertNextCell (int type, vtkIdType npts, vtkIdType *ptIds, vtkIdType nfaces, vtkIdType *faces)
 
int GetCellType (vtkIdType cellId)
 
vtkUnsignedCharArrayGetCellTypesArray ()
 
vtkIdTypeArrayGetCellLocationsArray ()
 
void Squeeze ()
 
void Initialize ()
 
int GetMaxCellSize ()
 
void BuildLinks ()
 
vtkCellLinksGetCellLinks ()
 
virtual void GetCellPoints (vtkIdType cellId, vtkIdType &npts, vtkIdType *&pts)
 
void GetFaceStream (vtkIdType cellId, vtkIdList *ptIds)
 
void GetFaceStream (vtkIdType cellId, vtkIdType &nfaces, vtkIdType *&ptIds)
 
vtkCellArrayGetCells ()
 
void ReplaceCell (vtkIdType cellId, int npts, vtkIdType *pts)
 
vtkIdType InsertNextLinkedCell (int type, int npts, vtkIdType *pts)
 
void RemoveReferenceToCell (vtkIdType ptId, vtkIdType cellId)
 
void AddReferenceToCell (vtkIdType ptId, vtkIdType cellId)
 
void ResizeCellList (vtkIdType ptId, int size)
 
virtual int GetGhostLevel ()
 
unsigned long GetActualMemorySize ()
 
void GetIdsOfCellsOfType (int type, vtkIdTypeArray *array)
 
int IsHomogeneous ()
 
void RemoveGhostCells ()
 
vtkIdTypeGetFaces (vtkIdType cellId)
 
int InitializeFacesRepresentation (vtkIdType numPrevCells)
 
void Reset ()
 
virtual void CopyStructure (vtkDataSet *ds)
 
vtkIdType GetNumberOfCells ()
 
virtual vtkCellGetCell (vtkIdType cellId)
 
virtual void GetCell (vtkIdType cellId, vtkGenericCell *cell)
 
virtual void GetCellBounds (vtkIdType cellId, double bounds[6])
 
virtual void GetCellPoints (vtkIdType cellId, vtkIdList *ptIds)
 
void GetPointCells (vtkIdType ptId, vtkIdList *cellIds)
 
vtkCellIteratorNewCellIterator ()
 
void SetCells (int type, vtkCellArray *cells)
 
void SetCells (int *types, vtkCellArray *cells)
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells)
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells, vtkIdTypeArray *faceLocations, vtkIdTypeArray *faces)
 
virtual void GetCellNeighbors (vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
 
virtual int GetPiece ()
 
virtual int GetNumberOfPieces ()
 
vtkIdTypeArrayGetFaces ()
 
vtkIdTypeArrayGetFaceLocations ()
 
void Reset ()
 
virtual void CopyStructure (vtkDataSet *ds)
 
vtkIdType GetNumberOfCells ()
 
virtual vtkCellGetCell (vtkIdType cellId)
 
virtual void GetCell (vtkIdType cellId, vtkGenericCell *cell)
 
virtual void GetCellBounds (vtkIdType cellId, double bounds[6])
 
virtual void GetCellPoints (vtkIdType cellId, vtkIdList *ptIds)
 
void GetPointCells (vtkIdType ptId, vtkIdList *cellIds)
 
vtkCellIteratorNewCellIterator ()
 
void SetCells (int type, vtkCellArray *cells)
 
void SetCells (int *types, vtkCellArray *cells)
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells)
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells, vtkIdTypeArray *faceLocations, vtkIdTypeArray *faces)
 
virtual void GetCellNeighbors (vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
 
virtual int GetPiece ()
 
virtual int GetNumberOfPieces ()
 
vtkIdTypeArrayGetFaces ()
 
vtkIdTypeArrayGetFaceLocations ()
 
- Public Member Functions inherited from vtkUnstructuredGridBase
 vtkAbstractTypeMacro (vtkUnstructuredGridBase, vtkPointSet) void PrintSelf(ostream &os
 
int GetDataObjectType ()
 
- Public Member Functions inherited from vtkPointSet
vtkPointSetNewInstance () const
 
doubleGetPoint (vtkIdType ptId)
 
vtkMTimeType GetMTime ()
 
void ComputeBounds ()
 
vtkIdType GetNumberOfPoints ()
 
void GetPoint (vtkIdType ptId, double x[3])
 
virtual vtkIdType FindPoint (double x[3])
 
vtkIdType FindPoint (double x, double y, double z)
 
virtual vtkIdType FindCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual vtkIdType FindCell (double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual void SetPoints (vtkPoints *)
 
virtual vtkPointsGetPoints ()
 
void Register (vtkObjectBase *o) VTK_OVERRIDE
 
void UnRegister (vtkObjectBase *o) VTK_OVERRIDE
 
vtkIdType GetNumberOfPoints ()
 
void GetPoint (vtkIdType ptId, double x[3])
 
virtual vtkIdType FindPoint (double x[3])
 
vtkIdType FindPoint (double x, double y, double z)
 
virtual vtkIdType FindCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual vtkIdType FindCell (double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual void SetPoints (vtkPoints *)
 
virtual vtkPointsGetPoints ()
 
void Register (vtkObjectBase *o) VTK_OVERRIDE
 
void UnRegister (vtkObjectBase *o) VTK_OVERRIDE
 
- Public Member Functions inherited from vtkDataSet
vtkDataSetNewInstance () const
 
virtual void CopyAttributes (vtkDataSet *ds)
 
virtual void GetCellTypes (vtkCellTypes *types)
 
vtkCellDataGetCellData ()
 
vtkPointDataGetPointData ()
 
doubleGetBounds ()
 
void GetBounds (double bounds[6])
 
doubleGetCenter ()
 
void GetCenter (double center[3])
 
double GetLength ()
 
virtual void GetScalarRange (double range[2])
 
doubleGetScalarRange ()
 
int CheckAttributes ()
 
virtual vtkFieldDataGetAttributesAsFieldData (int type)
 
virtual vtkIdType GetNumberOfElements (int type)
 
bool HasAnyGhostCells ()
 
bool HasAnyGhostPoints ()
 
vtkUnsignedCharArrayGetPointGhostArray ()
 
void UpdatePointGhostArrayCache ()
 
vtkUnsignedCharArrayAllocatePointGhostArray ()
 
vtkUnsignedCharArrayGetCellGhostArray ()
 
void UpdateCellGhostArrayCache ()
 
vtkUnsignedCharArrayAllocateCellGhostArray ()
 
vtkIdType FindPoint (double x, double y, double z)
 
virtual vtkCellFindAndGetCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual void GenerateGhostArray (int zeroExt[6])
 
virtual void GenerateGhostArray (int zeroExt[6], bool cellOnly)
 
virtual bool HasAnyBlankCells ()
 
virtual bool HasAnyBlankPoints ()
 
vtkIdType FindPoint (double x, double y, double z)
 
virtual vtkCellFindAndGetCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 
virtual void GenerateGhostArray (int zeroExt[6])
 
virtual void GenerateGhostArray (int zeroExt[6], bool cellOnly)
 
virtual bool HasAnyBlankCells ()
 
virtual bool HasAnyBlankPoints ()
 
- Public Member Functions inherited from vtkDataObject
vtkDataObjectNewInstance () const
 
void ReleaseData ()
 
vtkMTimeType GetUpdateTime ()
 
virtual void CopyInformationToPipeline (vtkInformation *vtkNotUsed(info))
 
void DataHasBeenGenerated ()
 
virtual void PrepareForNewData ()
 
virtual int GetExtentType ()
 
virtual void Crop (const int *updateExtent)
 
virtual vtkDataSetAttributesGetAttributes (int type)
 
virtual int GetAttributeTypeForArray (vtkAbstractArray *arr)
 
virtual vtkInformationGetInformation ()
 
virtual void SetInformation (vtkInformation *)
 
virtual int GetDataReleased ()
 
virtual void SetFieldData (vtkFieldData *)
 
virtual vtkFieldDataGetFieldData ()
 
virtual void CopyInformationFromPipeline (vtkInformation *vtkNotUsed(info))
 
virtual vtkInformationGetInformation ()
 
virtual void SetInformation (vtkInformation *)
 
virtual int GetDataReleased ()
 
void GlobalReleaseDataFlagOn ()
 
void GlobalReleaseDataFlagOff ()
 
virtual void SetFieldData (vtkFieldData *)
 
virtual vtkFieldDataGetFieldData ()
 
virtual void CopyInformationFromPipeline (vtkInformation *vtkNotUsed(info))
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 
virtual void DebugOff ()
 
bool GetDebug ()
 
void SetDebug (bool debugFlag)
 
virtual void Modified ()
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
int HasObserver (unsigned long event)
 
int HasObserver (const char *event)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 
vtkCommandGetCommand (unsigned long tag)
 
void RemoveObserver (vtkCommand *)
 
void RemoveObservers (unsigned long event, vtkCommand *)
 
void RemoveObservers (const char *event, vtkCommand *)
 
int HasObserver (unsigned long event, vtkCommand *)
 
int HasObserver (const char *event, vtkCommand *)
 
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
int InvokeEvent (unsigned long event, void *callData)
 
int InvokeEvent (const char *event, void *callData)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 
vtkCommandGetCommand (unsigned long tag)
 
void RemoveObserver (vtkCommand *)
 
void RemoveObservers (unsigned long event, vtkCommand *)
 
void RemoveObservers (const char *event, vtkCommand *)
 
int HasObserver (unsigned long event, vtkCommand *)
 
int HasObserver (const char *event, vtkCommand *)
 
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
int InvokeEvent (unsigned long event, void *callData)
 
int InvokeEvent (const char *event, void *callData)
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 
virtual void Delete ()
 
virtual void FastDelete ()
 
void Print (ostream &os)
 
void SetReferenceCount (int)
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 
int GetReferenceCount ()
 
void PrintRevisions (ostream &)
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 
int GetReferenceCount ()
 
void PrintRevisions (ostream &)
 

Static Public Member Functions

static vtkboneFiniteElementModelNew ()
 
static int IsTypeOf (const char *type)
 
static vtkboneFiniteElementModelSafeDownCast (vtkObject *o)
 
static const char * GetElementTypeAsString (int arg)
 
- Static Public Member Functions inherited from vtkUnstructuredGrid
static vtkUnstructuredGridNew ()
 
static int IsTypeOf (const char *type)
 
static vtkUnstructuredGridSafeDownCast (vtkObjectBase *o)
 
static void DecomposeAPolyhedronCell (vtkIdType *polyhedronCellStream, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
static vtkUnstructuredGridGetData (vtkInformation *info)
 
static vtkUnstructuredGridGetData (vtkInformationVector *v, int i=0)
 
static void DecomposeAPolyhedronCell (vtkCellArray *polyhedronCellArray, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
static void DecomposeAPolyhedronCell (vtkIdType nCellFaces, vtkIdType *inFaceStream, vtkIdType &nCellpts, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
static void ConvertFaceStreamPointIds (vtkIdList *faceStream, vtkIdType *idMap)
 
static void ConvertFaceStreamPointIds (vtkIdType nfaces, vtkIdType *faceStream, vtkIdType *idMap)
 
static vtkUnstructuredGridGetData (vtkInformation *info)
 
static vtkUnstructuredGridGetData (vtkInformationVector *v, int i=0)
 
static void DecomposeAPolyhedronCell (vtkCellArray *polyhedronCellArray, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
static void DecomposeAPolyhedronCell (vtkIdType nCellFaces, vtkIdType *inFaceStream, vtkIdType &nCellpts, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
static void ConvertFaceStreamPointIds (vtkIdList *faceStream, vtkIdType *idMap)
 
static void ConvertFaceStreamPointIds (vtkIdType nfaces, vtkIdType *faceStream, vtkIdType *idMap)
 
- Static Public Member Functions inherited from vtkUnstructuredGridBase
static vtkUnstructuredGridBaseGetData (vtkInformation *info)
 
static vtkUnstructuredGridBaseGetData (vtkInformationVector *v, int i=0)
 
static vtkUnstructuredGridBaseGetData (vtkInformation *info)
 
static vtkUnstructuredGridBaseGetData (vtkInformationVector *v, int i=0)
 
- Static Public Member Functions inherited from vtkPointSet
static int IsTypeOf (const char *type)
 
static vtkPointSetSafeDownCast (vtkObjectBase *o)
 
static vtkPointSetGetData (vtkInformation *info)
 
static vtkPointSetGetData (vtkInformationVector *v, int i=0)
 
static vtkPointSetGetData (vtkInformation *info)
 
static vtkPointSetGetData (vtkInformationVector *v, int i=0)
 
- Static Public Member Functions inherited from vtkDataSet
static int IsTypeOf (const char *type)
 
static vtkDataSetSafeDownCast (vtkObjectBase *o)
 
static vtkDataSetGetData (vtkInformation *info)
 
static vtkDataSetGetData (vtkInformationVector *v, int i=0)
 
static vtkDataSetGetData (vtkInformation *info)
 
static vtkDataSetGetData (vtkInformationVector *v, int i=0)
 
- Static Public Member Functions inherited from vtkDataObject
static vtkDataObjectNew ()
 
static int IsTypeOf (const char *type)
 
static vtkDataObjectSafeDownCast (vtkObjectBase *o)
 
static const char * GetAssociationTypeAsString (int associationType)
 
static int GetAssociationTypeFromString (const char *associationType)
 
static vtkInformationStringKeyDATA_TYPE_NAME ()
 
static vtkInformationDataObjectKeyDATA_OBJECT ()
 
static vtkInformationIntegerKeyDATA_EXTENT_TYPE ()
 
static vtkInformationIntegerPointerKeyDATA_EXTENT ()
 
static vtkInformationIntegerVectorKeyALL_PIECES_EXTENT ()
 
static vtkInformationIntegerKeyDATA_PIECE_NUMBER ()
 
static vtkInformationIntegerKeyDATA_NUMBER_OF_PIECES ()
 
static vtkInformationIntegerKeyDATA_NUMBER_OF_GHOST_LEVELS ()
 
static vtkInformationDoubleKeyDATA_TIME_STEP ()
 
static vtkInformationInformationVectorKeyPOINT_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyCELL_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyVERTEX_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyEDGE_DATA_VECTOR ()
 
static vtkInformationIntegerKeyFIELD_ARRAY_TYPE ()
 
static vtkInformationIntegerKeyFIELD_ASSOCIATION ()
 
static vtkInformationIntegerKeyFIELD_ATTRIBUTE_TYPE ()
 
static vtkInformationIntegerKeyFIELD_ACTIVE_ATTRIBUTE ()
 
static vtkInformationIntegerKeyFIELD_NUMBER_OF_COMPONENTS ()
 
static vtkInformationIntegerKeyFIELD_NUMBER_OF_TUPLES ()
 
static vtkInformationIntegerKeyFIELD_OPERATION ()
 
static vtkInformationDoubleVectorKeyFIELD_RANGE ()
 
static vtkInformationIntegerVectorKeyPIECE_EXTENT ()
 
static vtkInformationStringKeyFIELD_NAME ()
 
static vtkInformationDoubleVectorKeyORIGIN ()
 
static vtkInformationDoubleVectorKeySPACING ()
 
static vtkInformationDoubleVectorKeyBOUNDING_BOX ()
 
static vtkInformationDataObjectKeySIL ()
 
static vtkInformationGetActiveFieldInformation (vtkInformation *info, int fieldAssociation, int attributeType)
 
static vtkInformationGetNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static void RemoveNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static vtkInformationSetActiveAttribute (vtkInformation *info, int fieldAssociation, const char *attributeName, int attributeType)
 
static void SetActiveAttributeInfo (vtkInformation *info, int fieldAssociation, int attributeType, const char *name, int arrayType, int numComponents, int numTuples)
 
static void SetPointDataActiveScalarInfo (vtkInformation *info, int arrayType, int numComponents)
 
static vtkDataObjectGetData (vtkInformation *info)
 
static vtkDataObjectGetData (vtkInformationVector *v, int i=0)
 
static void SetGlobalReleaseDataFlag (int val)
 
static int GetGlobalReleaseDataFlag ()
 
static vtkInformationGetActiveFieldInformation (vtkInformation *info, int fieldAssociation, int attributeType)
 
static vtkInformationGetNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static void RemoveNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static vtkInformationSetActiveAttribute (vtkInformation *info, int fieldAssociation, const char *attributeName, int attributeType)
 
static void SetActiveAttributeInfo (vtkInformation *info, int fieldAssociation, int attributeType, const char *name, int arrayType, int numComponents, int numTuples)
 
static void SetPointDataActiveScalarInfo (vtkInformation *info, int arrayType, int numComponents)
 
static vtkDataObjectGetData (vtkInformation *info)
 
static vtkDataObjectGetData (vtkInformationVector *v, int i=0)
 
- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 
static void BreakOnError ()
 
static void SetGlobalWarningDisplay (int val)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static int GetGlobalWarningDisplay ()
 
static void SetGlobalWarningDisplay (int val)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static int GetGlobalWarningDisplay ()
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 
static vtkObjectBaseNew ()
 
static vtkObjectBaseNew ()
 

Protected Member Functions

 vtkboneFiniteElementModel ()
 
 ~vtkboneFiniteElementModel ()
 
- Protected Member Functions inherited from vtkUnstructuredGrid
virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkUnstructuredGrid ()
 
 ~vtkUnstructuredGrid ()
 
- Protected Member Functions inherited from vtkUnstructuredGridBase
 vtkUnstructuredGridBase ()
 
 ~vtkUnstructuredGridBase ()
 
- Protected Member Functions inherited from vtkPointSet
 vtkPointSet ()
 
 ~vtkPointSet ()
 
void ReportReferences (vtkGarbageCollector *) VTK_OVERRIDE
 
- Protected Member Functions inherited from vtkDataSet
 vtkDataSet ()
 
 ~vtkDataSet ()
 
virtual void ComputeScalarRange ()
 
bool IsAnyBitSet (vtkUnsignedCharArray *a, int bitFlag)
 
- Protected Member Functions inherited from vtkDataObject
 vtkDataObject ()
 
 ~vtkDataObject ()
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
virtual ~vtkObject ()
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) VTK_OVERRIDE
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) VTK_OVERRIDE
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=NULL)
 
void InternalReleaseFocus ()
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=NULL)
 
void InternalReleaseFocus ()
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkDataArrayCollectionNodeSets
 
vtkDataArrayCollectionElementSets
 
vtkDataArrayCollectionGaussPointData
 
vtkboneConstraintCollectionConstraints
 
vtkboneMaterialTableMaterialTable
 
vtkboneConstraintConvergenceSet
 
char * Name
 
char * History
 
char * Log
 
- Protected Attributes inherited from vtkUnstructuredGrid
vtkVertexVertex
 
vtkPolyVertexPolyVertex
 
vtkLineLine
 
vtkPolyLinePolyLine
 
vtkTriangleTriangle
 
vtkTriangleStripTriangleStrip
 
vtkPixelPixel
 
vtkQuadQuad
 
vtkPolygonPolygon
 
vtkTetraTetra
 
vtkVoxelVoxel
 
vtkHexahedronHexahedron
 
vtkWedgeWedge
 
vtkPyramidPyramid
 
vtkPentagonalPrismPentagonalPrism
 
vtkHexagonalPrismHexagonalPrism
 
vtkQuadraticEdgeQuadraticEdge
 
vtkQuadraticTriangleQuadraticTriangle
 
vtkQuadraticQuadQuadraticQuad
 
vtkQuadraticPolygonQuadraticPolygon
 
vtkQuadraticTetraQuadraticTetra
 
vtkQuadraticHexahedronQuadraticHexahedron
 
vtkQuadraticWedgeQuadraticWedge
 
vtkQuadraticPyramidQuadraticPyramid
 
vtkQuadraticLinearQuadQuadraticLinearQuad
 
vtkBiQuadraticQuadBiQuadraticQuad
 
vtkTriQuadraticHexahedronTriQuadraticHexahedron
 
vtkQuadraticLinearWedgeQuadraticLinearWedge
 
vtkBiQuadraticQuadraticWedgeBiQuadraticQuadraticWedge
 
vtkBiQuadraticQuadraticHexahedronBiQuadraticQuadraticHexahedron
 
vtkBiQuadraticTriangleBiQuadraticTriangle
 
vtkCubicLineCubicLine
 
vtkConvexPointSetConvexPointSet
 
vtkPolyhedronPolyhedron
 
vtkEmptyCellEmptyCell
 
vtkCellArrayConnectivity
 
vtkCellLinksLinks
 
vtkUnsignedCharArrayTypes
 
vtkIdTypeArrayLocations
 
vtkIdTypeArrayFaces
 
vtkIdTypeArrayFaceLocations
 
- Protected Attributes inherited from vtkPointSet
vtkPointsPoints
 
vtkPointLocatorLocator
 
- Protected Attributes inherited from vtkDataSet
vtkCellDataCellData
 
vtkPointDataPointData
 
vtkCallbackCommandDataObserver
 
vtkTimeStamp ComputeTime
 
double Bounds [6]
 
double Center [3]
 
double ScalarRange [2]
 
vtkTimeStamp ScalarRangeComputeTime
 
vtkUnsignedCharArrayPointGhostArray
 
vtkUnsignedCharArrayCellGhostArray
 
bool PointGhostArrayCached
 
bool CellGhostArrayCached
 
- Protected Attributes inherited from vtkDataObject
vtkFieldDataFieldData
 
int DataReleased
 
vtkTimeStamp UpdateTime
 
vtkInformationInformation
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkUnstructuredGridBase
vtkIndent indent
 
- Public Attributes inherited from vtkDataSet
 DATA_OBJECT_FIELD
 
 POINT_DATA_FIELD
 
 CELL_DATA_FIELD
 
- Public Attributes inherited from vtkDataObject
 FIELD_ASSOCIATION_POINTS
 
 FIELD_ASSOCIATION_CELLS
 
 FIELD_ASSOCIATION_NONE
 
 FIELD_ASSOCIATION_POINTS_THEN_CELLS
 
 FIELD_ASSOCIATION_VERTICES
 
 FIELD_ASSOCIATION_EDGES
 
 FIELD_ASSOCIATION_ROWS
 
 NUMBER_OF_ASSOCIATIONS
 
 POINT
 
 CELL
 
 FIELD
 
 POINT_THEN_CELL
 
 VERTEX
 
 EDGE
 
 ROW
 
 NUMBER_OF_ATTRIBUTE_TYPES
 
 FIELD_OPERATION_PRESERVED
 
 FIELD_OPERATION_REINTERPOLATED
 
 FIELD_OPERATION_MODIFIED
 
 FIELD_OPERATION_REMOVED
 

Detailed Description

data model for finite element meshes

vtkboneFiniteElementModel represents a finite element model. vtkboneFiniteElementModel is a subclass of vtkUnstructuredGrid. The underlying vtkUnstructuredGrid is the complete set of nodes and elements of the finite element model and thus specifies the geometry.

Named subsets of nodes and elements may be defined.

Optionally contains constraints (such as boundary conditions and applied loads), stored as vtkboneConstraintCollection.

Optionally contains a material property table, stored as vtkboneMaterialTable.

This object is 0-indexed on all arrays. Where output is required to be 1-indexed, translation is performed in the appropriate writer object (e.g. vtkboneN88ModelWriter).

The topology is as specified by VTK. For convenience, the topologies of the commonly used VTK_VOXEL and VTK_HEXAHEDRON are shown below.

  VTK_VOXEL
          6---------7
         /|        /|
        / |       / |
       /  |      /  |
      4---------5   |
      |   |     |   |
      |   2-----|---3
      |  /      |  /
      |/        | /     z  y
      0---------1/      | /
                        |/
                        .--->x

  VTK_HEXAHEDRON
          7---------6
         /|        /|
        / |       / |
       /  |      /  |
      4---------5   |
      |   |     |   |
      |   3-----|---2
      |  /      |  /
      |/        | /     z  y
      0---------1/      | /
                        |/
                        .--->x
See also
vtkboneConstraint vtkboneSolverParameters vtkboneMaterialTable vtkboneFiniteElementModelGenerator

Definition at line 87 of file vtkboneFiniteElementModel.h.

Member Enumeration Documentation

Specifies the type of elements in the model.

Enumerator
N88_UNKNOWN 
N88_TETRAHEDRON 
N88_HEXAHEDRON 
N88_MIXED 
NUMBER_OF_ElementType 

Definition at line 191 of file vtkboneFiniteElementModel.h.

Constructor & Destructor Documentation

vtkboneFiniteElementModel::vtkboneFiniteElementModel ( )
protected
vtkboneFiniteElementModel::~vtkboneFiniteElementModel ( )
protected

Member Function Documentation

static vtkboneFiniteElementModel* vtkboneFiniteElementModel::New ( )
static
virtual const char* vtkboneFiniteElementModel::GetClassName ( )
virtual
static int vtkboneFiniteElementModel::IsTypeOf ( const char *  type)
static
virtual int vtkboneFiniteElementModel::IsA ( const char *  type)
virtual

Reimplemented from vtkUnstructuredGrid.

static vtkboneFiniteElementModel* vtkboneFiniteElementModel::SafeDownCast ( vtkObject o)
static
void vtkboneFiniteElementModel::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
virtual

Reimplemented from vtkUnstructuredGrid.

virtual void vtkboneFiniteElementModel::AddNodeSet ( vtkIdTypeArray nodeSet)
virtual

Add a node set. Note that the node set must have a name or an error will occur. Any existing node set with the same name will be replaced.

virtual void vtkboneFiniteElementModel::AddElementSet ( vtkIdTypeArray elementSet)
virtual

Add an element set. Note that the element set must have a name or an error will occur. Any existing element set with the same name will be replaced.

virtual vtkIdTypeArray* vtkboneFiniteElementModel::GetNodeSet ( const char *  nodeSetName)
virtual

Returns a pointer to the named node set. Returns NULL if no such node set exists.

virtual vtkIdTypeArray* vtkboneFiniteElementModel::GetElementSet ( const char *  elementSetName)
virtual

Returns a pointer to the named element set. Returns NULL if no such element set exists.

virtual int vtkboneFiniteElementModel::RemoveNodeSet ( const char *  nodeSetName)
virtual

Remove a node set. Returns VTK_ERROR if the named node set does not exist.

virtual int vtkboneFiniteElementModel::RemoveElementSet ( const char *  elementSetName)
virtual

Remove an element set. Returns VTK_ERROR if the named element set does not exist.

virtual int vtkboneFiniteElementModel::GetAssociatedElementsFromNodeSet ( const char *  nodeSetName,
vtkIdTypeArray ids 
)
virtual

Returns the cell Ids of all elements that contain at least one point identified by the specified node set. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkIdTypeArray* vtkboneFiniteElementModel::GetAssociatedElementsFromNodeSet ( const char *  nodeSetName)
virtual

Returns the cell Ids of all elements that contain at least one point identified by the specified node set. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::DataSetFromNodeSet ( const vtkIdTypeArray nodeSet,
vtkUnstructuredGrid data 
)
virtual

Return a vtkUnstructuredGrid corresponding to the specified NodeSet. The result is a dataset of vertices on the selected nodes. This is useful for rendering the NodeSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points in the node set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromNodeSet ( const vtkIdTypeArray nodeSet)
virtual

Return a vtkUnstructuredGrid corresponding to the specified NodeSet. The result is a dataset of vertices on the selected nodes. This is useful for rendering the NodeSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points in the node set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::DataSetFromNodeSet ( const char *  nodeSetName,
vtkUnstructuredGrid data 
)
virtual

Return a vtkUnstructuredGrid corresponding to the specified NodeSet. The result is a dataset of vertices on the selected nodes. This is useful for rendering the NodeSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points in the node set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromNodeSet ( const char *  nodeSetName)
virtual

Return a vtkUnstructuredGrid corresponding to the specified NodeSet. The result is a dataset of vertices on the selected nodes. This is useful for rendering the NodeSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points in the node set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no node set with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::DataSetFromElementSet ( const vtkIdTypeArray elementSet,
vtkUnstructuredGrid data 
)
virtual

Return a vtkUnstructuredGrid corresponding to the specified ElementSet. The result is a subset of the cells of the model. This is useful for rendering the ElementSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points or cells in the element set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromElementSet ( const vtkIdTypeArray elementSet)
virtual

Return a vtkUnstructuredGrid corresponding to the specified ElementSet. The result is a subset of the cells of the model. This is useful for rendering the ElementSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points or cells in the element set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::DataSetFromElementSet ( const char *  elementSetName,
vtkUnstructuredGrid data 
)
virtual

Return a vtkUnstructuredGrid corresponding to the specified ElementSet. The result is a subset of the cells of the model. This is useful for rendering the ElementSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points or cells in the element set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromElementSet ( const char *  elementSetName)
virtual

Return a vtkUnstructuredGrid corresponding to the specified ElementSet. The result is a subset of the cells of the model. This is useful for rendering the ElementSet. Note: The Points list of the returned object is NOT the same as the vtkboneFiniteElementModel from which it is derived. However, you can still trace back the original node (point) ids from the PointData attribute PedigreeIds array, providing the vtkboneFiniteElementModel has such an array (by default it does). WARNING: Duplicate points or cells in the element set will be merged. The first version returns VTK_ERROR on error. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual void vtkboneFiniteElementModel::SetConstraints ( vtkboneConstraintCollection )
virtual

Set/get the constraints.

virtual vtkboneConstraintCollection* vtkboneFiniteElementModel::GetConstraints ( )
virtual

Set/get the constraints.

virtual void vtkboneFiniteElementModel::SetMaterialTable ( vtkboneMaterialTable )
virtual

Set/get the material table.

virtual vtkboneMaterialTable* vtkboneFiniteElementModel::GetMaterialTable ( )
virtual

Set/get the material table.

virtual void vtkboneFiniteElementModel::SetConvergenceSet ( vtkboneConstraint )
virtual

Set/get the convergence set. This is optional.

virtual vtkboneConstraint* vtkboneFiniteElementModel::GetConvergenceSet ( )
virtual

Set/get the convergence set. This is optional.

static const char* vtkboneFiniteElementModel::GetElementTypeAsString ( int  arg)
static

Return a string describing the value of ElementType_t

virtual int vtkboneFiniteElementModel::GetElementType ( )
virtual

Determine the correct element type. This call can be computationally expensive, so avoid calling it several times if possible.

virtual void vtkboneFiniteElementModel::GetAllCellPoints ( vtkIdTypeArray allCellPoints)
virtual

A convenience function that returns all cell points in one array. This is usually only useful if all cells are of the same type. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual vtkIdTypeArray* vtkboneFiniteElementModel::GetAllCellPoints ( )
virtual

A convenience function that returns all cell points in one array. This is usually only useful if all cells are of the same type. The second version returns NULL if no element set with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::ApplyBoundaryCondition ( vtkIdTypeArray nodeIds,
vtkDataArray senses,
vtkDataArray displacements,
const char *  arg_constraintName 
)
virtual

Generates a displacement Constraint along the specified axis sense with the specified amount of displacement. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyBoundaryCondition ( vtkIdTypeArray nodeIds,
int  sense,
double  displacement,
const char *  arg_constraintName 
)
virtual

Generates a displacement Constraint along the specified axis sense with the specified amount of displacement. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyBoundaryCondition ( vtkIdType  nodeId,
int  sense,
double  displacement,
const char *  arg_constraintName 
)
virtual

Generates a displacement Constraint along the specified axis sense with the specified amount of displacement. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyBoundaryCondition ( const char *  nodeSetName,
int  sense,
double  displacement,
const char *  arg_constraintName 
)
virtual

Generates a displacement Constraint along the specified axis sense with the specified amount of displacement. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::FixNodes ( vtkIdTypeArray ids,
const char *  arg_constraintName 
)
virtual

Generates a fixed Constraint, with displacement 0 applied to all senses of the specified nodes. This is a convenience function, as the same result can be obtained by calling ApplyBoundaryCondition for each axis, with displacement zero. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::FixNodes ( vtkIdType  id,
const char *  arg_constraintName 
)
virtual

Generates a fixed Constraint, with displacement 0 applied to all senses of the specified nodes. This is a convenience function, as the same result can be obtained by calling ApplyBoundaryCondition for each axis, with displacement zero. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::FixNodes ( const char *  selectionName,
const char *  arg_constraintName 
)
virtual

Generates a fixed Constraint, with displacement 0 applied to all senses of the specified nodes. This is a convenience function, as the same result can be obtained by calling ApplyBoundaryCondition for each axis, with displacement zero. The variations of this method allow the nodes to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyLoad ( vtkIdTypeArray elementIds,
vtkDataArray distributions,
vtkDataArray senses,
vtkDataArray forces,
const char *  arg_constraintName 
)
virtual

Generates a force Constraint along the specified axis sense with the specified amount of force. IMPORTANT: For the variants of this method that take a scalar value for the force, this value is the total force, which will be evenly distributed input set. The variations of this method allow the elements to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The distribution of the force on the elements must be specified with the parameter distribution; which must take a value of the enum Distribution_t defined in vtkboneConstraint. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyLoad ( vtkIdTypeArray elementIds,
int  distribution,
vtkDataArray senses,
vtkDataArray forces,
const char *  arg_constraintName 
)
virtual

Generates a force Constraint along the specified axis sense with the specified amount of force. IMPORTANT: For the variants of this method that take a scalar value for the force, this value is the total force, which will be evenly distributed input set. The variations of this method allow the elements to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The distribution of the force on the elements must be specified with the parameter distribution; which must take a value of the enum Distribution_t defined in vtkboneConstraint. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyLoad ( vtkIdTypeArray elementIds,
int  distribution,
int  sense,
double  total_force,
const char *  arg_constraintName 
)
virtual

Generates a force Constraint along the specified axis sense with the specified amount of force. IMPORTANT: For the variants of this method that take a scalar value for the force, this value is the total force, which will be evenly distributed input set. The variations of this method allow the elements to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The distribution of the force on the elements must be specified with the parameter distribution; which must take a value of the enum Distribution_t defined in vtkboneConstraint. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyLoad ( vtkIdType  elementId,
int  distribution,
int  sense,
double  force,
const char *  arg_constraintName 
)
virtual

Generates a force Constraint along the specified axis sense with the specified amount of force. IMPORTANT: For the variants of this method that take a scalar value for the force, this value is the total force, which will be evenly distributed input set. The variations of this method allow the elements to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The distribution of the force on the elements must be specified with the parameter distribution; which must take a value of the enum Distribution_t defined in vtkboneConstraint. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ApplyLoad ( const char *  elementSetName,
int  distribution,
int  sense,
double  total_force,
const char *  arg_constraintName 
)
virtual

Generates a force Constraint along the specified axis sense with the specified amount of force. IMPORTANT: For the variants of this method that take a scalar value for the force, this value is the total force, which will be evenly distributed input set. The variations of this method allow the elements to be specified as a vtkIdTypeArray list of node Ids, as a single node Id, or as a named existing node set. The distribution of the force on the elements must be specified with the parameter distribution; which must take a value of the enum Distribution_t defined in vtkboneConstraint. The Constraint is added to this object's constraint collection with the specified name. If an existing constraint has the same name, then the new constraints are added to it.

virtual int vtkboneFiniteElementModel::ConvergenceSetFromConstraint ( vtkboneConstraint constraint)
virtual

Given a constraint, generate the complementary convergence set from it. For example, for a constraint on displacements, a convergence set of force will be generated.

virtual int vtkboneFiniteElementModel::DataSetFromConstraint ( vtkboneConstraint arg_constraint,
vtkUnstructuredGrid data 
)
virtual

Return a Dataset consisting of the nodes or cells associated with the constraint. If the constraint is a cell constraint, you will get that subset of cells. If the constraint is a node constraint, you will get a dataset of vertices on the constrained nodes. The first version returns VTK_ERROR on error. The second version returns NULL if no constraint with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromConstraint ( vtkboneConstraint arg_constraint)
virtual

Return a Dataset consisting of the nodes or cells associated with the constraint. If the constraint is a cell constraint, you will get that subset of cells. If the constraint is a node constraint, you will get a dataset of vertices on the constrained nodes. The first version returns VTK_ERROR on error. The second version returns NULL if no constraint with the specified name exists. You are responsible for deleting the returned object.

virtual int vtkboneFiniteElementModel::DataSetFromConstraint ( const char *  arg_constraintName,
vtkUnstructuredGrid data 
)
virtual

Return a Dataset consisting of the nodes or cells associated with the constraint. If the constraint is a cell constraint, you will get that subset of cells. If the constraint is a node constraint, you will get a dataset of vertices on the constrained nodes. The first version returns VTK_ERROR on error. The second version returns NULL if no constraint with the specified name exists. You are responsible for deleting the returned object.

virtual vtkUnstructuredGrid* vtkboneFiniteElementModel::DataSetFromConstraint ( const char *  arg_constraintName)
virtual

Return a Dataset consisting of the nodes or cells associated with the constraint. If the constraint is a cell constraint, you will get that subset of cells. If the constraint is a node constraint, you will get a dataset of vertices on the constrained nodes. The first version returns VTK_ERROR on error. The second version returns NULL if no constraint with the specified name exists. You are responsible for deleting the returned object.

virtual void vtkboneFiniteElementModel::SetName ( const char *  )
virtual

Set/get the name of the model.

virtual char* vtkboneFiniteElementModel::GetName ( )
virtual

Set/get the name of the model.

virtual void vtkboneFiniteElementModel::SetNodeSets ( vtkDataArrayCollection )
virtual

Set/get the node sets.

virtual vtkDataArrayCollection* vtkboneFiniteElementModel::GetNodeSets ( )
virtual

Set/get the node sets.

virtual void vtkboneFiniteElementModel::SetElementSets ( vtkDataArrayCollection )
virtual

Set/get the element sets.

virtual vtkDataArrayCollection* vtkboneFiniteElementModel::GetElementSets ( )
virtual

Set/get the element sets.

virtual void vtkboneFiniteElementModel::SetGaussPointData ( vtkDataArrayCollection )
virtual

Set/get the gauss point data. Each set of gauss point data should be stored as a named vtkFloatArray. The number of components of the array should be set to the number of values per Gauss point (e.g. 1 for scalar data, 6 for stress/strain data). The number of tuples is then the product of the number of elements (i.e. Cells) times the number of gauss points per element.

virtual vtkDataArrayCollection* vtkboneFiniteElementModel::GetGaussPointData ( )
virtual

Set/get the gauss point data. Each set of gauss point data should be stored as a named vtkFloatArray. The number of components of the array should be set to the number of values per Gauss point (e.g. 1 for scalar data, 6 for stress/strain data). The number of tuples is then the product of the number of elements (i.e. Cells) times the number of gauss points per element.

virtual void vtkboneFiniteElementModel::SetHistory ( const char *  )
virtual

Set/Get the History. The history should be a list of one-liners, each line starting with the data and followed by the user and/or program and/or object that is creating or modifying the model. If more information is required, add it to the Log field.

virtual char* vtkboneFiniteElementModel::GetHistory ( )
virtual

Set/Get the History. The history should be a list of one-liners, each line starting with the data and followed by the user and/or program and/or object that is creating or modifying the model. If more information is required, add it to the Log field.

void vtkboneFiniteElementModel::AppendHistory ( const char *  s)

Append to the History. The date will automatically be added.

virtual void vtkboneFiniteElementModel::SetLog ( const char *  )
virtual

Detailed information on how the model was created and/or modified.

virtual char* vtkboneFiniteElementModel::GetLog ( )
virtual

Detailed information on how the model was created and/or modified.

void vtkboneFiniteElementModel::AppendLog ( const char *  s)

Append to the Log. The date will automatically be added.

virtual void vtkboneFiniteElementModel::ShallowCopy ( vtkDataObject src)
virtual

Shallow and Deep copy.

Reimplemented from vtkUnstructuredGrid.

virtual void vtkboneFiniteElementModel::DeepCopy ( vtkDataObject src)
virtual

Shallow and Deep copy.

Reimplemented from vtkUnstructuredGrid.

Member Data Documentation

vtkDataArrayCollection* vtkboneFiniteElementModel::NodeSets
protected

Definition at line 385 of file vtkboneFiniteElementModel.h.

vtkDataArrayCollection* vtkboneFiniteElementModel::ElementSets
protected

Definition at line 386 of file vtkboneFiniteElementModel.h.

vtkDataArrayCollection* vtkboneFiniteElementModel::GaussPointData
protected

Definition at line 387 of file vtkboneFiniteElementModel.h.

vtkboneConstraintCollection* vtkboneFiniteElementModel::Constraints
protected

Definition at line 389 of file vtkboneFiniteElementModel.h.

vtkboneMaterialTable* vtkboneFiniteElementModel::MaterialTable
protected

Definition at line 390 of file vtkboneFiniteElementModel.h.

vtkboneConstraint* vtkboneFiniteElementModel::ConvergenceSet
protected

Definition at line 391 of file vtkboneFiniteElementModel.h.

char* vtkboneFiniteElementModel::Name
protected

Definition at line 393 of file vtkboneFiniteElementModel.h.

char* vtkboneFiniteElementModel::History
protected

Definition at line 394 of file vtkboneFiniteElementModel.h.

char* vtkboneFiniteElementModel::Log
protected

Definition at line 395 of file vtkboneFiniteElementModel.h.


The documentation for this class was generated from the following file: