VTK
vtkQuadratureSchemeDefinition.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadratureSchemeDefinition.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
35 #ifndef vtkQuadratureSchemeDefinition_h
36 #define vtkQuadratureSchemeDefinition_h
37 
38 #include "vtkCommonDataModelModule.h" // For export macro
39 #include "vtkObject.h"
40 
43 class vtkXMLDataElement;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkQuadratureSchemeDefinition : public vtkObject
46 {
47 public:
48  // vtk stuff
50  void PrintSelf(ostream& os, vtkIndent indent) override;
52  static vtkInformationStringKey* QUADRATURE_OFFSET_ARRAY_NAME();
53 
59 
63  int DeepCopy(const vtkQuadratureSchemeDefinition *other);
64 
69  int SaveState(vtkXMLDataElement *e);
73  int RestoreState(vtkXMLDataElement *e);
74 
79  void Clear();
80 
84  void Initialize(int cellType,
85  int numberOfNodes,
86  int numberOfQuadraturePoints,
87  double *shapeFunctionWeights);
91  void Initialize(int cellType,
92  int numberOfNodes,
93  int numberOfQuadraturePoints,
94  double *shapeFunctionWeights,
95  double *quadratureWeights);
96 
100  int GetCellType() const { return this->CellType; }
104  int GetQuadratureKey() const { return this->QuadratureKey; }
108  int GetNumberOfNodes() const { return this->NumberOfNodes; }
112  int GetNumberOfQuadraturePoints() const { return this->NumberOfQuadraturePoints; }
118  const double *GetShapeFunctionWeights() const { return this->ShapeFunctionWeights; }
120 
124  const double *GetShapeFunctionWeights(int quadraturePointId) const
125  {
126  int idx=quadraturePointId*this->NumberOfNodes;
127  return this->ShapeFunctionWeights+idx;
128  }
130 
133  const double *GetQuadratureWeights() const { return this->QuadratureWeights; }
134 
135 protected:
137  ~vtkQuadratureSchemeDefinition() override;
138 private:
143  void ReleaseResources();
148  int SecureResources();
153  void SetShapeFunctionWeights(const double *W);
158  void SetQuadratureWeights(const double *W);
159 
160  //
162  void operator=(const vtkQuadratureSchemeDefinition &) = delete;
163  friend ostream &operator<<(ostream &s, const vtkQuadratureSchemeDefinition &d);
164  friend istream &operator>>(istream &s, vtkQuadratureSchemeDefinition &d);
165  //
166  int CellType;
167  int QuadratureKey;
168  int NumberOfNodes;
169  int NumberOfQuadraturePoints;
170  double *ShapeFunctionWeights;
171  double *QuadratureWeights;
172 };
173 
174 #endif
175 
abstract base class for most VTK objects
Definition: vtkObject.h:59
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
Represents an XML element and those nested inside.
int GetQuadratureKey() const
Access to an alternative key.
int GetNumberOfQuadraturePoints() const
Get the number of quadrature points associated with the scheme.
An Elemental data type that holds a definition of a numerical quadrature scheme.
const double * GetShapeFunctionWeights() const
Get the array of shape function weights.
Key for string values in vtkInformation.
int GetCellType() const
Access the VTK cell type id.
int GetNumberOfNodes() const
Get the number of nodes associated with the interpolation.
a simple class to control print indentation
Definition: vtkIndent.h:39
const double * GetQuadratureWeights() const
Access to the quadrature weights.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
friend VTKCOMMONCORE_EXPORT ostream & operator<<(ostream &os, vtkObjectBase &o)
const double * GetShapeFunctionWeights(int quadraturePointId) const
Get the array of shape function weights associated with a single quadrature point.