VTK  9.5.2
vtkQuadratureSchemeDefinition.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
23#ifndef vtkQuadratureSchemeDefinition_h
24#define vtkQuadratureSchemeDefinition_h
25
26#include "vtkCellType.h" // For VTK_EMPTY_CELL
27#include "vtkCommonDataModelModule.h" // For export macro
28#include "vtkObject.h"
29
30VTK_ABI_NAMESPACE_BEGIN
34
35class VTKCOMMONDATAMODEL_EXPORT vtkQuadratureSchemeDefinition : public vtkObject
36{
37public:
38 // vtk stuff
40 void PrintSelf(ostream& os, vtkIndent indent) override;
43
49
54
60
65
70 int cellType, int numberOfNodes, int numberOfQuadraturePoints, double* shapeFunctionWeights);
71
75 void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints,
76 double* shapeFunctionWeights, double* quadratureWeights);
77
81 void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints,
82 const double* shapeFunctionWeights, const double* quadratureWeights, int dim,
83 const double* shapeFunctionDerivativeWeights);
84
88 int GetCellType() const { return this->CellType; }
89
93 int GetQuadratureKey() const { return this->QuadratureKey; }
94
98 int GetNumberOfNodes() const { return this->NumberOfNodes; }
99
103 int GetNumberOfQuadraturePoints() const { return this->NumberOfQuadraturePoints; }
104
108 vtkGetMacro(Dimension, int);
109
115 const double* GetShapeFunctionWeights() const { return this->ShapeFunctionWeights; }
116
118
122 const double* GetShapeFunctionWeights(int quadraturePointId) const
123 {
124 int idx = quadraturePointId * this->NumberOfNodes;
125 return this->ShapeFunctionWeights + idx;
126 }
127
132 const double* GetShapeFunctionDerivativeWeights(int quadraturePointId) const
133 {
134 int idx = quadraturePointId * this->NumberOfNodes * this->Dimension;
135 return this->ShapeFunctionDerivativeWeights + idx;
136 }
138
142 const double* GetQuadratureWeights() const { return this->QuadratureWeights; }
143
144protected:
147
148private:
153 void ReleaseResources();
154
159 int SecureResources();
160
165 void SetShapeFunctionWeights(const double* weights);
166
171 void SetQuadratureWeights(const double* weights);
172
177 void SetShapeFunctionDerivativeWeights(const double* weights);
178
179 //
181 void operator=(const vtkQuadratureSchemeDefinition&) = delete;
182 friend ostream& operator<<(ostream& s, const vtkQuadratureSchemeDefinition& d);
183 friend istream& operator>>(istream& s, vtkQuadratureSchemeDefinition& d);
184 //
185 int CellType = VTK_EMPTY_CELL;
186 int QuadratureKey = -1;
187 int NumberOfNodes = 0;
188 int NumberOfQuadraturePoints = 0;
189 int Dimension = 0;
190 double* ShapeFunctionWeights = nullptr;
191 double* QuadratureWeights = nullptr;
192 double* ShapeFunctionDerivativeWeights = nullptr;
193};
194
195VTK_ABI_NAMESPACE_END
196#endif
a simple class to control print indentation
Definition vtkIndent.h:108
Key for string values in vtkInformation.
abstract base class for most VTK objects
Definition vtkObject.h:162
An Elemental data type that holds a definition of a numerical quadrature scheme.
friend ostream & operator<<(ostream &s, const vtkQuadratureSchemeDefinition &d)
void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints, double *shapeFunctionWeights)
Initialize the object allocating resources as needed.
~vtkQuadratureSchemeDefinition() override
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetQuadratureKey() const
Access to an alternative key.
void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints, const double *shapeFunctionWeights, const double *quadratureWeights, int dim, const double *shapeFunctionDerivativeWeights)
Initialize the object allocating resources as needed.
int GetNumberOfQuadraturePoints() const
Get the number of quadrature points associated with the scheme.
static vtkInformationStringKey * QUADRATURE_OFFSET_ARRAY_NAME()
static vtkInformationQuadratureSchemeDefinitionVectorKey * DICTIONARY()
void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints, double *shapeFunctionWeights, double *quadratureWeights)
Initialize the object allocating resources as needed.
int RestoreState(vtkXMLDataElement *root)
Restore the object from an XML representation.
static vtkQuadratureSchemeDefinition * New()
New object in an unusable state.
const double * GetQuadratureWeights() const
Access to the quadrature weights.
const double * GetShapeFunctionWeights(int quadraturePointId) const
Get the array of shape function weights associated with a single quadrature point.
int SaveState(vtkXMLDataElement *root)
Put the object into an XML representation.
int DeepCopy(const vtkQuadratureSchemeDefinition *other)
Deep copy.
const double * GetShapeFunctionDerivativeWeights(int quadraturePointId) const
Get the array of shape function derivative weights associated with a single quadrature point.
friend istream & operator>>(istream &s, vtkQuadratureSchemeDefinition &d)
const double * GetShapeFunctionWeights() const
Get the array of shape function weights.
int GetCellType() const
Access the VTK cell type id.
int GetNumberOfNodes() const
Get the number of nodes associated with the interpolation.
Represents an XML element and those nested inside.
@ VTK_EMPTY_CELL
Definition vtkCellType.h:37