VTK  9.1.0
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
44
45class VTKCOMMONDATAMODEL_EXPORT vtkQuadratureSchemeDefinition : public vtkObject
46{
47public:
48 // vtk stuff
50 void PrintSelf(ostream& os, vtkIndent indent) override;
53
59
64
74
79 void Clear();
80
85 int cellType, int numberOfNodes, int numberOfQuadraturePoints, double* shapeFunctionWeights);
89 void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints,
90 double* shapeFunctionWeights, double* quadratureWeights);
91
95 int GetCellType() const { return this->CellType; }
99 int GetQuadratureKey() const { return this->QuadratureKey; }
103 int GetNumberOfNodes() const { return this->NumberOfNodes; }
107 int GetNumberOfQuadraturePoints() const { return this->NumberOfQuadraturePoints; }
113 const double* GetShapeFunctionWeights() const { return this->ShapeFunctionWeights; }
115
119 const double* GetShapeFunctionWeights(int quadraturePointId) const
120 {
121 int idx = quadraturePointId * this->NumberOfNodes;
122 return this->ShapeFunctionWeights + idx;
123 }
125
128 const double* GetQuadratureWeights() const { return this->QuadratureWeights; }
129
130protected:
133
134private:
139 void ReleaseResources();
144 int SecureResources();
149 void SetShapeFunctionWeights(const double* W);
154 void SetQuadratureWeights(const double* W);
155
156 //
158 void operator=(const vtkQuadratureSchemeDefinition&) = delete;
159 friend ostream& operator<<(ostream& s, const vtkQuadratureSchemeDefinition& d);
160 friend istream& operator>>(istream& s, vtkQuadratureSchemeDefinition& d);
161 //
162 int CellType;
163 int QuadratureKey;
164 int NumberOfNodes;
165 int NumberOfQuadraturePoints;
166 double* ShapeFunctionWeights;
167 double* QuadratureWeights;
168};
169
170#endif
a simple class to control print indentation
Definition: vtkIndent.h:113
Key for string values in vtkInformation.
abstract base class for most VTK objects
Definition: vtkObject.h:73
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.
void Clear()
Release all allocated resources and set the object to an uninitialized state.
~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.
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 unsuable 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.
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.