vtkbone
vtkboneFiniteElementModel.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Copyright 2010-2016, Numerics88 Solutions Ltd.
4 http://www.numerics88.com/
5
6 Copyright (c) Eric Nodwell and Steven K. Boyd
7 See Copyright.txt for details.
8
9 This software is distributed WITHOUT ANY WARRANTY; without even
10 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11 PURPOSE. See the above copyright notice for more information.
12=========================================================================*/
13
72#ifndef __vtkboneFiniteElementModel_h
73#define __vtkboneFiniteElementModel_h
74
75#include "vtkUnstructuredGrid.h"
76#include "vtkboneWin32Header.h"
77#include "vtkSmartPointer.h"
78
79// Forward declarations
80class vtkIdTypeArray;
81class vtkCharArray;
86
88{
89public:
92 void PrintSelf(ostream& os, vtkIndent indent) override;
93
97 virtual void AddNodeSet (vtkIdTypeArray* nodeSet);
98
102 virtual void AddElementSet (vtkIdTypeArray* elementSet);
103
106 virtual vtkIdTypeArray* GetNodeSet (const char* nodeSetName);
107
110 virtual vtkIdTypeArray* GetElementSet (const char* elementSetName);
111
114 virtual int RemoveNodeSet (const char* nodeSetName);
115
118 virtual int RemoveElementSet (const char* elementSetName);
119
121
126 virtual int GetAssociatedElementsFromNodeSet (const char *nodeSetName,
127 vtkIdTypeArray* ids);
128 virtual vtkIdTypeArray* GetAssociatedElementsFromNodeSet (const char *nodeSetName);
130
132
143 virtual int DataSetFromNodeSet (const vtkIdTypeArray *nodeSet,
144 vtkUnstructuredGrid* data);
146 virtual int DataSetFromNodeSet (const char *nodeSetName,
147 vtkUnstructuredGrid* data);
148 virtual vtkUnstructuredGrid* DataSetFromNodeSet (const char *nodeSetName);
150
152
163 virtual int DataSetFromElementSet(const vtkIdTypeArray *elementSet,
164 vtkUnstructuredGrid* data);
166 virtual int DataSetFromElementSet(const char *elementSetName,
167 vtkUnstructuredGrid* data);
168 virtual vtkUnstructuredGrid* DataSetFromElementSet(const char *elementSetName);
170
172
174 vtkGetObjectMacro(Constraints, vtkboneConstraintCollection);
176
178
180 vtkGetObjectMacro(MaterialTable, vtkboneMaterialTable);
182
184
186 vtkGetObjectMacro(ConvergenceSet, vtkboneConstraint);
188
190
196 NUMBER_OF_ElementType};
198
199 static const char* GetElementTypeAsString (int arg);
200
203 virtual int GetElementType();
204
206
210 virtual void GetAllCellPoints (vtkIdTypeArray* allCellPoints);
213
215
223 vtkIdTypeArray* nodeIds,
224 vtkDataArray* senses,
225 vtkDataArray* displacements,
226 const char* arg_constraintName);
228 vtkIdTypeArray* nodeIds,
229 int sense,
230 double displacement,
231 const char* arg_constraintName);
233 vtkIdType nodeId,
234 int sense,
235 double displacement,
236 const char* arg_constraintName);
238 const char* nodeSetName,
239 int sense,
240 double displacement,
241 const char* arg_constraintName);
243
245
254 virtual int FixNodes(vtkIdTypeArray* ids, const char* arg_constraintName);
255 virtual int FixNodes(vtkIdType id, const char* arg_constraintName);
256 virtual int FixNodes(const char* selectionName, const char* arg_constraintName);
258
260
272 virtual int ApplyLoad(
273 vtkIdTypeArray* elementIds,
274 vtkDataArray* distributions,
275 vtkDataArray* senses,
276 vtkDataArray* forces,
277 const char* arg_constraintName);
278 virtual int ApplyLoad(
279 vtkIdTypeArray* elementIds,
280 int distribution,
281 vtkDataArray* senses,
282 vtkDataArray* forces,
283 const char* arg_constraintName);
284 virtual int ApplyLoad(
285 vtkIdTypeArray* elementIds,
286 int distribution,
287 int sense,
288 double total_force,
289 const char* arg_constraintName);
290 virtual int ApplyLoad(
291 vtkIdType elementId,
292 int distribution,
293 int sense,
294 double force,
295 const char* arg_constraintName);
296 virtual int ApplyLoad(
297 const char* elementSetName,
298 int distribution,
299 int sense,
300 double total_force,
301 const char* arg_constraintName);
303
305
309 virtual int ConvergenceSetFromConstraint(const char* constraintName);
311
313
320 virtual int DataSetFromConstraint(vtkboneConstraint* arg_constraint,
321 vtkUnstructuredGrid* data);
323 virtual int DataSetFromConstraint(const char *arg_constraintName,
324 vtkUnstructuredGrid* data);
325 virtual vtkUnstructuredGrid* DataSetFromConstraint(const char *arg_constraintName);
327
329
330 vtkSetStringMacro(Name);
331 vtkGetStringMacro(Name);
333
335
337 vtkGetObjectMacro(NodeSets, vtkDataArrayCollection);
339
341
343 vtkGetObjectMacro(ElementSets, vtkDataArrayCollection);
345
347
354 vtkGetObjectMacro(GaussPointData, vtkDataArrayCollection);
356
358
362 vtkSetStringMacro(History);
363 vtkGetStringMacro(History);
365
367 void AppendHistory(const char* s);
368
370
371 vtkSetStringMacro(Log);
372 vtkGetStringMacro(Log);
374
376 void AppendLog(const char* s);
377
379
380 virtual void ShallowCopy(vtkDataObject *src) override;
381 virtual void DeepCopy(vtkDataObject *src) override;
383
384protected:
387
391
395
396 char* Name;
397 char* History;
398 char* Log;
399
400private:
401 vtkboneFiniteElementModel(const vtkboneFiniteElementModel&); // Not implemented.
402 void operator=(const vtkboneFiniteElementModel&); // Not implemented.
403};
404
405#endif
406
void operator=(const vtkObjectBase &)
maintain an unordered list of dataarray objects
a constraint for a finite element mesh
data model for finite element meshes
virtual int ApplyBoundaryCondition(vtkIdTypeArray *nodeIds, vtkDataArray *senses, vtkDataArray *displacements, const char *arg_constraintName)
virtual void ShallowCopy(vtkDataObject *src) override
virtual int FixNodes(const char *selectionName, const char *arg_constraintName)
virtual int DataSetFromElementSet(const char *elementSetName, vtkUnstructuredGrid *data)
static const char * GetElementTypeAsString(int arg)
virtual vtkUnstructuredGrid * DataSetFromConstraint(vtkboneConstraint *arg_constraint)
virtual int ApplyLoad(vtkIdType elementId, int distribution, int sense, double force, const char *arg_constraintName)
virtual int RemoveElementSet(const char *elementSetName)
virtual int ApplyLoad(vtkIdTypeArray *elementIds, vtkDataArray *distributions, vtkDataArray *senses, vtkDataArray *forces, const char *arg_constraintName)
virtual void DeepCopy(vtkDataObject *src) override
virtual void AddNodeSet(vtkIdTypeArray *nodeSet)
virtual int FixNodes(vtkIdTypeArray *ids, const char *arg_constraintName)
vtkDataArrayCollection * NodeSets
virtual int DataSetFromNodeSet(const vtkIdTypeArray *nodeSet, vtkUnstructuredGrid *data)
virtual void SetConstraints(vtkboneConstraintCollection *)
virtual int DataSetFromElementSet(const vtkIdTypeArray *elementSet, vtkUnstructuredGrid *data)
virtual vtkIdTypeArray * GetAssociatedElementsFromNodeSet(const char *nodeSetName)
virtual int RemoveNodeSet(const char *nodeSetName)
virtual void SetMaterialTable(vtkboneMaterialTable *)
virtual int ApplyLoad(vtkIdTypeArray *elementIds, int distribution, int sense, double total_force, const char *arg_constraintName)
void PrintSelf(ostream &os, vtkIndent indent) override
vtkboneConstraintCollection * Constraints
virtual vtkIdTypeArray * GetElementSet(const char *elementSetName)
virtual int ApplyBoundaryCondition(vtkIdType nodeId, int sense, double displacement, const char *arg_constraintName)
virtual int FixNodes(vtkIdType id, const char *arg_constraintName)
virtual void SetGaussPointData(vtkDataArrayCollection *)
vtkDataArrayCollection * GaussPointData
virtual int ApplyBoundaryCondition(vtkIdTypeArray *nodeIds, int sense, double displacement, const char *arg_constraintName)
virtual vtkUnstructuredGrid * DataSetFromConstraint(const char *arg_constraintName)
virtual int DataSetFromConstraint(const char *arg_constraintName, vtkUnstructuredGrid *data)
virtual int DataSetFromNodeSet(const char *nodeSetName, vtkUnstructuredGrid *data)
virtual vtkIdTypeArray * GetNodeSet(const char *nodeSetName)
vtkDataArrayCollection * ElementSets
static vtkboneFiniteElementModel * New()
virtual int DataSetFromConstraint(vtkboneConstraint *arg_constraint, vtkUnstructuredGrid *data)
void AppendLog(const char *s)
virtual int GetAssociatedElementsFromNodeSet(const char *nodeSetName, vtkIdTypeArray *ids)
virtual void SetNodeSets(vtkDataArrayCollection *)
virtual void SetElementSets(vtkDataArrayCollection *)
void AppendHistory(const char *s)
virtual void SetConvergenceSet(vtkboneConstraint *)
virtual int GetElementType()
virtual void GetAllCellPoints(vtkIdTypeArray *allCellPoints)
virtual int ConvergenceSetFromConstraint(const char *constraintName)
virtual int ApplyLoad(const char *elementSetName, int distribution, int sense, double total_force, const char *arg_constraintName)
virtual vtkUnstructuredGrid * DataSetFromNodeSet(const char *nodeSetName)
virtual vtkUnstructuredGrid * DataSetFromElementSet(const vtkIdTypeArray *elementSet)
virtual vtkIdTypeArray * GetAllCellPoints()
virtual int ApplyBoundaryCondition(const char *nodeSetName, int sense, double displacement, const char *arg_constraintName)
virtual vtkUnstructuredGrid * DataSetFromElementSet(const char *elementSetName)
virtual vtkUnstructuredGrid * DataSetFromNodeSet(const vtkIdTypeArray *nodeSet)
virtual void AddElementSet(vtkIdTypeArray *elementSet)
virtual int ApplyLoad(vtkIdTypeArray *elementIds, int distribution, vtkDataArray *senses, vtkDataArray *forces, const char *arg_constraintName)
virtual int ConvergenceSetFromConstraint(vtkboneConstraint *constraint)
Material Table finite element mesh.
int vtkIdType