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;
53 
59 
64 
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:
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 
vtkInformationStringKey
Key for string values in vtkInformation.
Definition: vtkInformationStringKey.h:37
vtkQuadratureSchemeDefinition::DICTIONARY
static vtkInformationQuadratureSchemeDefinitionVectorKey * DICTIONARY()
vtkQuadratureSchemeDefinition::operator>>
friend istream & operator>>(istream &s, vtkQuadratureSchemeDefinition &d)
vtkQuadratureSchemeDefinition::GetQuadratureWeights
const double * GetQuadratureWeights() const
Access to the quadrature weights.
Definition: vtkQuadratureSchemeDefinition.h:133
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:60
vtkQuadratureSchemeDefinition::DeepCopy
int DeepCopy(const vtkQuadratureSchemeDefinition *other)
Deep copy.
vtkQuadratureSchemeDefinition::operator<<
friend ostream & operator<<(ostream &s, const vtkQuadratureSchemeDefinition &d)
vtkQuadratureSchemeDefinition
Definition: vtkQuadratureSchemeDefinition.h:46
vtkQuadratureSchemeDefinition::New
static vtkQuadratureSchemeDefinition * New()
New object in an unsuable state.
vtkQuadratureSchemeDefinition::GetNumberOfNodes
int GetNumberOfNodes() const
Get the number of nodes associated with the interpolation.
Definition: vtkQuadratureSchemeDefinition.h:108
vtkQuadratureSchemeDefinition::~vtkQuadratureSchemeDefinition
~vtkQuadratureSchemeDefinition() override
vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME
static vtkInformationStringKey * QUADRATURE_OFFSET_ARRAY_NAME()
vtkQuadratureSchemeDefinition::GetCellType
int GetCellType() const
Access the VTK cell type id.
Definition: vtkQuadratureSchemeDefinition.h:100
vtkQuadratureSchemeDefinition::GetNumberOfQuadraturePoints
int GetNumberOfQuadraturePoints() const
Get the number of quadrature points associated with the scheme.
Definition: vtkQuadratureSchemeDefinition.h:112
vtkQuadratureSchemeDefinition::GetQuadratureKey
int GetQuadratureKey() const
Access to an alternative key.
Definition: vtkQuadratureSchemeDefinition.h:104
vtkQuadratureSchemeDefinition::GetShapeFunctionWeights
const double * GetShapeFunctionWeights() const
Get the array of shape function weights.
Definition: vtkQuadratureSchemeDefinition.h:118
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
vtkXMLDataElement
Represents an XML element and those nested inside.
Definition: vtkXMLDataElement.h:37
vtkQuadratureSchemeDefinition::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkQuadratureSchemeDefinition::vtkQuadratureSchemeDefinition
vtkQuadratureSchemeDefinition()
vtkQuadratureSchemeDefinition::Initialize
void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints, double *shapeFunctionWeights)
Initialize the object allocating resources as needed.
vtkObject.h
vtkQuadratureSchemeDefinition::Clear
void Clear()
Release all allocated resources and set the object to an uninitialized state.
vtkQuadratureSchemeDefinition::SaveState
int SaveState(vtkXMLDataElement *e)
Put the object into an XML representation.
vtkQuadratureSchemeDefinition::GetShapeFunctionWeights
const double * GetShapeFunctionWeights(int quadraturePointId) const
Get the array of shape function weights associated with a single quadrature point.
Definition: vtkQuadratureSchemeDefinition.h:124
vtkInformationQuadratureSchemeDefinitionVectorKey
Key for vtkQuadratureSchemeDefinition vector values.
Definition: vtkInformationQuadratureSchemeDefinitionVectorKey.h:37
vtkQuadratureSchemeDefinition::RestoreState
int RestoreState(vtkXMLDataElement *e)
Restore the object from an XML representation.
vtkQuadratureSchemeDefinition::Initialize
void Initialize(int cellType, int numberOfNodes, int numberOfQuadraturePoints, double *shapeFunctionWeights, double *quadratureWeights)
Initialize the object allocating resources as needed.