VTK  9.5.20250803
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;
41
43
49
55
60
66
71
76 int cellType, int numberOfNodes, int numberOfQuadraturePoints, double* shapeFunctionWeights);
77
81 void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints,
82 double* shapeFunctionWeights, double* quadratureWeights);
83
87 void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints,
88 const double* shapeFunctionWeights, const double* quadratureWeights, int dim,
89 const double* shapeFunctionDerivativeWeights);
90
94 int GetCellType() const { return this->CellType; }
95
99 int GetQuadratureKey() const { return this->QuadratureKey; }
100
104 int GetNumberOfNodes() const { return this->NumberOfNodes; }
105
109 int GetNumberOfQuadraturePoints() const { return this->NumberOfQuadraturePoints; }
110
114 vtkGetMacro(Dimension, int);
115
121 const double* GetShapeFunctionWeights() const { return this->ShapeFunctionWeights; }
122
124
128 const double* GetShapeFunctionWeights(int quadraturePointId) const
129 {
130 int idx = quadraturePointId * this->NumberOfNodes;
131 return this->ShapeFunctionWeights + idx;
132 }
133
138 const double* GetShapeFunctionDerivativeWeights(int quadraturePointId) const
139 {
140 int idx = quadraturePointId * this->NumberOfNodes * this->Dimension;
141 return this->ShapeFunctionDerivativeWeights + idx;
142 }
144
148 const double* GetQuadratureWeights() const { return this->QuadratureWeights; }
149
150protected:
153
154private:
159 void ReleaseResources();
160
165 int SecureResources();
166
171 void SetShapeFunctionWeights(const double* weights);
172
177 void SetQuadratureWeights(const double* weights);
178
183 void SetShapeFunctionDerivativeWeights(const double* weights);
184
185 //
187 void operator=(const vtkQuadratureSchemeDefinition&) = delete;
188 friend ostream& operator<<(ostream& s, const vtkQuadratureSchemeDefinition& d);
189 friend istream& operator>>(istream& s, vtkQuadratureSchemeDefinition& d);
190 //
191 int CellType = VTK_EMPTY_CELL;
192 int QuadratureKey = -1;
193 int NumberOfNodes = 0;
194 int NumberOfQuadraturePoints = 0;
195 int Dimension = 0;
196 double* ShapeFunctionWeights = nullptr;
197 double* QuadratureWeights = nullptr;
198 double* ShapeFunctionDerivativeWeights = nullptr;
199};
200
201VTK_ABI_NAMESPACE_END
202#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.
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.
static vtkInformationStringKey * QUADRATURE_OFFSET_ARRAY_NAME()
static vtkInformationQuadratureSchemeDefinitionVectorKey * DICTIONARY()
@ VTK_EMPTY_CELL
Definition vtkCellType.h:37