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
80 class vtkIdTypeArray;
81 class vtkCharArray;
83 class vtkboneConstraint;
86 
88 {
89 public:
92  void PrintSelf(ostream& os, vtkIndent indent);
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,
145  virtual vtkUnstructuredGrid* DataSetFromNodeSet (const vtkIdTypeArray *nodeSet);
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);
165  virtual vtkUnstructuredGrid* DataSetFromElementSet(const vtkIdTypeArray *elementSet);
166  virtual int DataSetFromElementSet(const char *elementSetName,
167  vtkUnstructuredGrid* data);
168  virtual vtkUnstructuredGrid* DataSetFromElementSet(const char *elementSetName);
170 
172 
173  virtual void SetConstraints(vtkboneConstraintCollection *);
174  vtkGetObjectMacro(Constraints, vtkboneConstraintCollection);
176 
178 
179  virtual void SetMaterialTable(vtkboneMaterialTable *);
180  vtkGetObjectMacro(MaterialTable, vtkboneMaterialTable);
182 
184 
185  virtual void SetConvergenceSet(vtkboneConstraint *);
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);
211  virtual vtkIdTypeArray* GetAllCellPoints ();
213 
215 
222  virtual int ApplyBoundaryCondition(
223  vtkIdTypeArray* nodeIds,
224  vtkDataArray* senses,
225  vtkDataArray* displacements,
226  const char* arg_constraintName);
227  virtual int ApplyBoundaryCondition(
228  vtkIdTypeArray* nodeIds,
229  int sense,
230  double displacement,
231  const char* arg_constraintName);
232  virtual int ApplyBoundaryCondition(
233  vtkIdType nodeId,
234  int sense,
235  double displacement,
236  const char* arg_constraintName);
237  virtual int ApplyBoundaryCondition(
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 
307  virtual int ConvergenceSetFromConstraint(vtkboneConstraint* constraint);
308 
310 
317  virtual int DataSetFromConstraint(vtkboneConstraint* arg_constraint,
318  vtkUnstructuredGrid* data);
319  virtual vtkUnstructuredGrid* DataSetFromConstraint(vtkboneConstraint* arg_constraint);
320  virtual int DataSetFromConstraint(const char *arg_constraintName,
321  vtkUnstructuredGrid* data);
322  virtual vtkUnstructuredGrid* DataSetFromConstraint(const char *arg_constraintName);
324 
326 
327  vtkSetStringMacro(Name);
328  vtkGetStringMacro(Name);
330 
332 
333  virtual void SetNodeSets (vtkDataArrayCollection*);
334  vtkGetObjectMacro(NodeSets, vtkDataArrayCollection);
336 
338 
339  virtual void SetElementSets (vtkDataArrayCollection*);
340  vtkGetObjectMacro(ElementSets, vtkDataArrayCollection);
342 
344 
350  virtual void SetGaussPointData(vtkDataArrayCollection *);
351  vtkGetObjectMacro(GaussPointData, vtkDataArrayCollection);
353 
355 
359  vtkSetStringMacro(History);
360  vtkGetStringMacro(History);
362 
364  void AppendHistory(const char* s);
365 
367 
368  vtkSetStringMacro(Log);
369  vtkGetStringMacro(Log);
371 
373  void AppendLog(const char* s);
374 
376 
377  virtual void ShallowCopy(vtkDataObject *src);
378  virtual void DeepCopy(vtkDataObject *src);
380 
381 protected:
384 
388 
392 
393  char* Name;
394  char* History;
395  char* Log;
396 
397 private:
398  vtkboneFiniteElementModel(const vtkboneFiniteElementModel&); // Not implemented.
399  void operator=(const vtkboneFiniteElementModel&); // Not implemented.
400 };
401 
402 #endif
403 
data model for finite element meshes
data
#define VTKBONE_EXPORT
vtkDataArrayCollection * NodeSets
virtual void DeepCopy(vtkDataObject *src)
void PrintSelf(ostream &os, vtkIndent indent) VTK_OVERRIDE
int vtkIdType
virtual void ShallowCopy(vtkDataObject *src)
a constraint for a finite element mesh
Material Table finite element mesh.
static vtkUnstructuredGrid * New()
vtkboneMaterialTable * MaterialTable
vtkboneConstraintCollection * Constraints
vtkDataArrayCollection * GaussPointData
vtkDataArrayCollection * ElementSets
maintain an unordered list of dataarray objects
void operator=(const vtkObjectBase &)