VTK  9.2.6
vtkDataSetEdgeSubdivisionCriterion.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkDataSetEdgeSubdivisionCriterion.h
5 Language: C++
6
7 Copyright 2003 Sandia Corporation.
8 Under the terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 license for use of this work by or on behalf of the
10 U.S. Government. Redistribution and use in source and binary forms, with
11 or without modification, are permitted provided that this Notice and any
12 statement of authorship are reproduced on all copies.
13
14=========================================================================*/
15#ifndef vtkDataSetEdgeSubdivisionCriterion_h
16#define vtkDataSetEdgeSubdivisionCriterion_h
40#include "vtkFiltersCoreModule.h" // For export macro
41
42class vtkCell;
43class vtkDataSet;
44
46{
47public:
50 void PrintSelf(ostream& os, vtkIndent indent) override;
51
52 virtual void SetMesh(vtkDataSet*);
53 vtkDataSet* GetMesh();
54
55 const vtkDataSet* GetMesh() const;
56
57 virtual void SetCellId(vtkIdType cell);
58 vtkIdType GetCellId() const;
59
60 vtkIdType& GetCellId();
61
62 vtkCell* GetCell();
63
64 const vtkCell* GetCell() const;
65
66 bool EvaluateLocationAndFields(double* midpt, int field_start) override;
67
106 double* EvaluateFields(double* vertex, double* weights, int field_start);
107
109
114 void EvaluatePointDataField(double* result, double* weights, int field);
115 void EvaluateCellDataField(double* result, double* weights, int field);
117
119
123 vtkSetMacro(ChordError2, double);
124 vtkGetMacro(ChordError2, double);
126
128
134 virtual void SetFieldError2(int s, double err);
135 double GetFieldError2(int s) const;
137
142 virtual void ResetFieldError2();
143
145
150 vtkGetMacro(ActiveFieldCriteria, int);
151
152// With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
153#if !VTK_USE_FUTURE_CONST
154 int GetActiveFieldCriteria() const { return this->ActiveFieldCriteria; }
155#endif
157
158protected:
161
165
167 double* FieldError2;
171
172private:
174 void operator=(const vtkDataSetEdgeSubdivisionCriterion&) = delete;
175};
176
178{
179 return this->CurrentCellId;
180}
182{
183 return this->CurrentCellId;
184}
185
187{
188 return this->CurrentMesh;
189}
191{
192 return this->CurrentMesh;
193}
194
196{
197 return this->CurrentCellData;
198}
200{
201 return this->CurrentCellData;
202}
203
204#endif // vtkDataSetEdgeSubdivisionCriterion_h
abstract class to specify cell behavior
Definition vtkCell.h:61
a subclass of vtkEdgeSubdivisionCriterion for vtkDataSet objects.
virtual void SetMesh(vtkDataSet *)
virtual void ResetFieldError2()
Tell the subdivider not to use any field values as subdivision criteria.
void EvaluateCellDataField(double *result, double *weights, int field)
Evaluate either a cell or nodal field.
virtual void SetFieldError2(int s, double err)
Get/Set the square of the allowable error magnitude for the scalar field s at any edge's midpoint.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
double * EvaluateFields(double *vertex, double *weights, int field_start)
Evaluate all of the fields that should be output with the given vertex and store them just past the p...
virtual void SetCellId(vtkIdType cell)
void EvaluatePointDataField(double *result, double *weights, int field)
Evaluate either a cell or nodal field.
bool EvaluateLocationAndFields(double *midpt, int field_start) override
You must implement this member function in a subclass.
static vtkDataSetEdgeSubdivisionCriterion * New()
int GetActiveFieldCriteria() const
Return a bitfield specifying which FieldError2 criteria are positive (i.e., actively used to decide e...
double GetFieldError2(int s) const
Get/Set the square of the allowable error magnitude for the scalar field s at any edge's midpoint.
abstract class to specify dataset behavior
Definition vtkDataSet.h:63
how to decide whether a linear approximation to nonlinear geometry or field should be subdivided
a simple class to control print indentation
Definition vtkIndent.h:40
int vtkIdType
Definition vtkType.h:332