gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
PViewData.h
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #ifndef PVIEW_DATA_H
7 #define PVIEW_DATA_H
8 
9 #include <string>
10 #include <vector>
11 #include <map>
12 #include <set>
13 #include "SBoundingBox3d.h"
14 #include "SPoint3KDTree.h"
15 
16 #define VAL_INF 1.e200
17 
18 class adaptiveData;
19 class GModel;
20 class GEntity;
21 class MElement;
22 class nameData;
23 class OctreePost;
24 template <class scalar> class fullMatrix;
25 
26 typedef std::map<int, std::vector<fullMatrix<double> *> > interpolationMatrices;
27 
28 // The abstract interface to post-processing view data.
29 class PViewData {
30 private:
31  // flag to mark that the data is 'dirty' and should not be displayed
32  bool _dirty;
33  // name of the view
34  std::string _name;
35  // name of the file the data was loaded from
36  std::string _fileName;
37  std::set<std::string> _fileNames;
38  // index of the view in the file
40  // octree for rapid search
42  // kdtree for rapid search of neighrest neighbor
46 
47 protected:
48  // adaptive visualization data
50  // interpolation matrices, indexed by the type of element
52  // global map of "named" interpolation matrices
53  static std::map<std::string, interpolationMatrices> _interpolationSchemes;
54  // string for the name of the interpolation scheme
56 
57 public:
58  PViewData();
59  virtual ~PViewData();
60 
61  // get/set the dirty ("not ready for display") flag
62  virtual bool getDirty() { return _dirty; }
63  virtual void setDirty(bool val) { _dirty = val; }
64 
65  // finalize the view data (compute min/max, etc.)
66  virtual bool finalize(bool computeMinMax = true,
67  const std::string &interpolationScheme = "");
68 
69  // get/set name
70  virtual std::string getName() { return _name; }
71  virtual void setName(const std::string &val) { _name = val; }
72 
73  // get/set (the main) filename containing the data
74  virtual std::string getFileName(int step = -1) { return _fileName; }
75  virtual void setFileName(const std::string &val)
76  {
77  _fileName = val;
78  _fileNames.insert(val);
79  }
80  virtual bool hasFileName(const std::string &val)
81  {
82  return _fileNames.find(val) != _fileNames.end();
83  }
84 
85  // get/set index of view data in file
86  virtual int getFileIndex() { return _fileIndex; }
87  virtual void setFileIndex(int val) { _fileIndex = val; }
88 
89  // get number of time steps in the data
90  virtual int getNumTimeSteps() = 0;
91  virtual int getFirstNonEmptyTimeStep(int start = 0) { return start; }
92 
93  // get the time value associated with the step-th time step
94  virtual double getTime(int step) { return 0.; }
95 
96  // get/set min/max for given step (global over all steps if step=-1)
97  virtual double getMin(int step = -1, bool onlyVisible = false,
98  int tensorRep = 0, int forceNumComponents = 0,
99  int componentMap[9] = nullptr) = 0;
100  virtual double getMax(int step = -1, bool onlyVisible = false,
101  int tensorRep = 0, int forceNumComponents = 0,
102  int componentMap[9] = nullptr) = 0;
103  virtual void setMin(double min) = 0;
104  virtual void setMax(double max) = 0;
105 
106  // get/set the bounding box
107  virtual SBoundingBox3d getBoundingBox(int step = -1) = 0;
108  virtual void setBoundingBox(SBoundingBox3d &box) = 0;
109 
110  // get the number of elements of a given type, for a given step
111  virtual int getNumScalars(int step = -1) { return 0; }
112  virtual int getNumVectors(int step = -1) { return 0; }
113  virtual int getNumTensors(int step = -1) { return 0; }
114  virtual int getNumPoints(int step = -1) { return 0; }
115  virtual int getNumLines(int step = -1) { return 0; }
116  virtual int getNumTriangles(int step = -1) { return 0; }
117  virtual int getNumQuadrangles(int step = -1) { return 0; }
118  virtual int getNumPolygons(int step = -1) { return 0; }
119  virtual int getNumTetrahedra(int step = -1) { return 0; }
120  virtual int getNumHexahedra(int step = -1) { return 0; }
121  virtual int getNumPrisms(int step = -1) { return 0; }
122  virtual int getNumPyramids(int step = -1) { return 0; }
123  virtual int getNumTrihedra(int step = -1) { return 0; }
124  virtual int getNumPolyhedra(int step = -1) { return 0; }
125 
126  // return the number of geometrical entities in the view
127  virtual int getNumEntities(int step = -1) { return 0; }
128 
129  // return the number of elements in the ent-th entity, or the total number of
130  // elements if ent < 0
131  virtual int getNumElements(int step = -1, int ent = -1) { return 0; }
132 
133  // return the geometrical dimension of the ele-th element in the ent-th entity
134  virtual int getDimension(int step, int ent, int ele) { return 0; }
135 
136  // return the number of nodes of the ele-th element in the ent-th entity
137  virtual int getNumNodes(int step, int ent, int ele) { return 0; }
138 
139  // get/set the coordinates and tag of the nod-th node from the ele-th element
140  // in the ent-th entity (if the node has a tag, getNode returns it)
141  virtual int getNode(int step, int ent, int ele, int nod, double &x, double &y,
142  double &z)
143  {
144  return 0;
145  }
146  virtual void setNode(int step, int ent, int ele, int nod, double x, double y,
147  double z);
148  virtual void tagNode(int step, int ent, int ele, int nod, int tag) {}
149 
150  // return the number of components available for the ele-th element in the
151  // ent-th entity
152  virtual int getNumComponents(int step, int ent, int ele) { return 0; }
153 
154  // return the number of values available for the ele-th element in the ent-th
155  // entity
156  virtual int getNumValues(int step, int ent, int ele) { return 0; }
157 
158  // get the idx'th value for the ele-th element in the ent-th entity
159  virtual void getValue(int step, int ent, int ele, int idx, double &val) {}
160 
161  // gets/set the comp-th component (at the step-th time step) associated with
162  // the node-th node from the ele-th element in the ent-th entity
163  virtual void getValue(int step, int ent, int ele, int nod, int comp,
164  double &val)
165  {
166  }
167  virtual void setValue(int step, int ent, int ele, int nod, int comp,
168  double val);
169 
170  // return a scalar value associated with the node-th node from the ele-th
171  // element in the ent-th entity: same as value for scalars, norm for vectors,
172  // Von-Mises (if tensorRep == 0), max eigenvalue (if tensorRep == 1) or min
173  // eigenvalue (if tensorRep == 2) for tensors
174  void getScalarValue(int step, int ent, int ele, int nod, double &val,
175  int tensorRep = 0, int forceNumComponents = 0,
176  int componentMap[9] = nullptr);
177 
178  // return the number of edges of the ele-th element in the ent-th
179  // entity
180  virtual int getNumEdges(int step, int ent, int ele) { return 0; }
181 
182  // return the type of the ele-th element in the ent-th entity
183  virtual int getType(int step, int ent, int ele) { return 0; }
184 
185  // return the number of 2D/3D strings in the view
186  virtual int getNumStrings2D() { return 0; }
187  virtual int getNumStrings3D() { return 0; }
188 
189  // return the i-th 2D/3D string in the view
190  virtual void getString2D(int i, int step, std::string &str, double &x,
191  double &y, double &style)
192  {
193  }
194  virtual void getString3D(int i, int step, std::string &str, double &x,
195  double &y, double &z, double &style)
196  {
197  }
198 
199  // change the orientation of the ele-th element
200  virtual void reverseElement(int step, int ent, int ele) {}
201 
202  // check if the view is empty
203  virtual bool empty();
204 
205  // check if we should skip the given entity/element
206  virtual bool skipEntity(int step, int ent) { return false; }
207  virtual bool skipElement(int step, int ent, int ele,
208  bool checkVisibility = false, int samplingRate = 1);
209 
210  // check if the data has the given step/partition/etc.
211  virtual bool hasTimeStep(int step)
212  {
213  return step >= 0 && step < getNumTimeSteps();
214  }
215  virtual bool hasPartition(int step, int part) { return false; }
216  virtual bool hasMultipleMeshes() { return false; }
217  virtual bool hasModel(GModel *model, int step = -1) { return false; }
218  virtual bool isNodeData() { return false; }
219 
220  // true if data is given at Gauss points (instead of vertices)
221  virtual bool useGaussPoints() { return false; }
222 
223  // initialize/destroy adaptive data
224  void initAdaptiveData(int step, int level, double tol);
225 
226  // Routines for
227  // - export of adapted views to pvtu file format for parallel visualization
228  // with paraview,
229  // - and/or generation of VTK data structure for ParaView plugin.
230  void initAdaptiveDataLight(int step, int level, double tol);
231  void saveAdaptedViewForVTK(const std::string &guifileName, int useDefaultName,
232  int step, int level, double tol, int npart,
233  bool isBinary);
234 
235  void destroyAdaptiveData();
236 
237  // return the adaptive data
239 
240  // set/get the interpolation matrices for elements with type "type"
241  // FIXME : too much overload :(
242  void setInterpolationMatrices(int type, const fullMatrix<double> &coefVal,
243  const fullMatrix<double> &expVal);
244  void setInterpolationMatrices(int type, const fullMatrix<double> &coefVal,
245  const fullMatrix<double> &expVal,
246  const fullMatrix<double> &coefGeo,
247  const fullMatrix<double> &expGeo);
248  int getInterpolationMatrices(int type, std::vector<fullMatrix<double> *> &p);
249  bool haveInterpolationMatrices(int type = 0);
250  void deleteInterpolationMatrices(int type = 0);
251 
252  // access to global interpolation schemes
253  static void removeInterpolationScheme(const std::string &name);
254  static void removeAllInterpolationSchemes();
255  static void addMatrixToInterpolationScheme(const std::string &name, int type,
256  fullMatrix<double> &mat);
257  static int getSizeInterpolationScheme();
259  void setInterpolationSchemeName(std::string name)
260  {
262  }
263 
264  // smooth the data in the view (makes it C0)
265  virtual void smooth();
266 
267  // combine time steps or elements from multiple datasets
268  virtual bool combineTime(nameData &nd);
269  virtual bool combineSpace(nameData &nd);
270 
271  // set simple X-Y data
272  virtual void setXY(std::vector<double> &x, std::vector<double> &y) {}
273 
274  // set simple pointwise XYZ data
275  virtual void setXYZV(std::vector<double> &x, std::vector<double> &y,
276  std::vector<double> &z, std::vector<double> &v)
277  {
278  }
279 
280  // ask to fill vertex arrays remotely
281  virtual bool isRemote() { return false; }
282  virtual int fillRemoteVertexArrays(std::string &options) { return 0; }
283 
284  // is the view a list-based dataset
285  virtual bool isListBased() { return false; }
286 
287  // get (approx) memry used by data in Mb
288  virtual double getMemoryInMb() { return 0; }
289 
290  // get GModel (if view supports it)
291  virtual GModel *getModel(int step);
292 
293  // get GEntity (if view supports it)
294  virtual GEntity *getEntity(int step, int entity);
295 
296  // get MElement (if view supports it)
297  virtual MElement *getElement(int step, int entity, int element);
298 
299  // find coordinates of closest node to point (xn, yn, zn); currently performs
300  // a simple linear search - we might want to use a kdtree instead
301  double findClosestNode(double &xn, double &yn, double &zn, int step);
302 
303  // search for the value of the View at point x, y, z. Values are interpolated
304  // using standard first order shape functions in the post element. If several
305  // time steps are present, they are all interpolated unless time step is set
306  // to a different value than -1.
307  bool searchScalar(double x, double y, double z, double *values, int step = -1,
308  double *size = nullptr, int qn = 0, double *qx = nullptr,
309  double *qy = nullptr, double *qz = nullptr,
310  bool grad = false, int dim = -1);
311  bool searchVector(double x, double y, double z, double *values, int step = -1,
312  double *size = nullptr, int qn = 0, double *qx = nullptr,
313  double *qy = nullptr, double *qz = nullptr,
314  bool grad = false, int dim = -1);
315  bool searchTensor(double x, double y, double z, double *values, int step = -1,
316  double *size = nullptr, int qn = 0, double *qx = nullptr,
317  double *qy = nullptr, double *qz = nullptr,
318  bool grad = false, int dim = -1);
319 
320  // same as above (when distance == 0), except that if no exact match is found:
321  // - if distance is > 0, return value at closest node if closer than distance
322  // - if distance is < 0, return value at closest node
323  // upon return, distance contains the distance of the match, or -1 if no
324  // match was found
325  bool searchScalarClosest(double x, double y, double z, double &distance,
326  double *values, int step = -1,
327  double *size = nullptr, int qn = 0, double *qx = nullptr,
328  double *qy = nullptr, double *qz = nullptr,
329  bool grad = false, int dim = -1);
330  bool searchVectorClosest(double x, double y, double z, double &distance,
331  double *values, int step = -1,
332  double *size = nullptr, int qn = 0, double *qx = nullptr,
333  double *qy = nullptr, double *qz = nullptr,
334  bool grad = false, int dim = -1);
335  bool searchTensorClosest(double x, double y, double z, double &distance,
336  double *values, int step = -1,
337  double *size = nullptr, int qn = 0, double *qx = nullptr,
338  double *qy = nullptr, double *qz = nullptr,
339  bool grad = false, int dim = -1);
340 
341  // I/O routines
342  virtual bool writeSTL(const std::string &fileName);
343  virtual bool writeTXT(const std::string &fileName);
344  virtual bool writePOS(const std::string &fileName, bool binary = false,
345  bool parsed = true, bool append = false);
346  virtual bool writeMSH(const std::string &fileName, double version = 2.2,
347  bool binary = false, bool saveMesh = true,
348  bool multipleView = false, int partitionNum = -1,
349  bool saveInterpolationMatrices = true,
350  bool forceNodeData = false,
351  bool forceElementData = false);
352  virtual bool writeMED(const std::string &fileName);
353  virtual bool toVector(std::vector<std::vector<double> > &vec);
354  virtual bool fromVector(const std::vector<std::vector<double> > &vec);
355  virtual void importLists(int N[24], std::vector<double> *V[24]);
356  virtual void getListPointers(int N[24], std::vector<double> *V[24]);
357  virtual void sendToServer(const std::string &name);
358 };
359 
360 #endif
PViewData::toVector
virtual bool toVector(std::vector< std::vector< double > > &vec)
Definition: PViewDataIO.cpp:217
PViewData::getNumEdges
virtual int getNumEdges(int step, int ent, int ele)
Definition: PViewData.h:180
PViewData::getNumPolygons
virtual int getNumPolygons(int step=-1)
Definition: PViewData.h:118
PViewData::getString3D
virtual void getString3D(int i, int step, std::string &str, double &x, double &y, double &z, double &style)
Definition: PViewData.h:194
PViewData::searchVector
bool searchVector(double x, double y, double z, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:362
PViewData::getNumLines
virtual int getNumLines(int step=-1)
Definition: PViewData.h:115
PViewData::writeSTL
virtual bool writeSTL(const std::string &fileName)
Definition: PViewDataIO.cpp:15
PViewData::getNumPyramids
virtual int getNumPyramids(int step=-1)
Definition: PViewData.h:122
PViewData::writeTXT
virtual bool writeTXT(const std::string &fileName)
Definition: PViewDataIO.cpp:74
PViewData::skipElement
virtual bool skipElement(int step, int ent, int ele, bool checkVisibility=false, int samplingRate=1)
Definition: PViewData.cpp:90
distance
double distance(MVertex *v1, MVertex *v2)
Definition: MVertex.h:245
PViewData::getNumStrings3D
virtual int getNumStrings3D()
Definition: PViewData.h:187
PViewData::combineSpace
virtual bool combineSpace(nameData &nd)
Definition: PViewData.cpp:242
PViewData::getNumStrings2D
virtual int getNumStrings2D()
Definition: PViewData.h:186
PViewData::getNumTimeSteps
virtual int getNumTimeSteps()=0
PViewData::setValue
virtual void setValue(int step, int ent, int ele, int nod, int comp, double val)
Definition: PViewData.cpp:130
PViewData::_name
std::string _name
Definition: PViewData.h:34
PViewData::isNodeData
virtual bool isNodeData()
Definition: PViewData.h:218
PViewData::getNode
virtual int getNode(int step, int ent, int ele, int nod, double &x, double &y, double &z)
Definition: PViewData.h:141
PViewData::_pc
SPoint3Cloud _pc
Definition: PViewData.h:43
PViewData::getInterpolationMatrices
int getInterpolationMatrices(int type, std::vector< fullMatrix< double > * > &p)
Definition: PViewData.cpp:176
box
Definition: gl2gif.cpp:311
SPoint3KDTree.h
PViewData::getValue
virtual void getValue(int step, int ent, int ele, int idx, double &val)
Definition: PViewData.h:159
PViewData::searchVectorClosest
bool searchVectorClosest(double x, double y, double z, double &distance, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:378
nameData
Definition: PView.h:169
OctreePost
Definition: BoundaryLayers.cpp:23
PViewData::getNumPoints
virtual int getNumPoints(int step=-1)
Definition: PViewData.h:114
PViewData::writePOS
virtual bool writePOS(const std::string &fileName, bool binary=false, bool parsed=true, bool append=false)
Definition: PViewDataIO.cpp:106
PViewData::_interpolationSchemeName
std::string _interpolationSchemeName
Definition: PViewData.h:55
PViewData::getNumValues
virtual int getNumValues(int step, int ent, int ele)
Definition: PViewData.h:156
PViewData::hasMultipleMeshes
virtual bool hasMultipleMeshes()
Definition: PViewData.h:216
PViewData::_octree
OctreePost * _octree
Definition: PViewData.h:41
PViewData::_interpolation
interpolationMatrices _interpolation
Definition: PViewData.h:51
PViewData::writeMED
virtual bool writeMED(const std::string &fileName)
Definition: PViewDataIO.cpp:211
PViewData::hasModel
virtual bool hasModel(GModel *model, int step=-1)
Definition: PViewData.h:217
PViewData::reverseElement
virtual void reverseElement(int step, int ent, int ele)
Definition: PViewData.h:200
PViewData::setNode
virtual void setNode(int step, int ent, int ele, int nod, double x, double y, double z)
Definition: PViewData.cpp:124
PViewData::getNumEntities
virtual int getNumEntities(int step=-1)
Definition: PViewData.h:127
PViewData::haveInterpolationMatrices
bool haveInterpolationMatrices(int type=0)
Definition: PViewData.cpp:186
GEntity
Definition: GEntity.h:31
PViewData::hasPartition
virtual bool hasPartition(int step, int part)
Definition: PViewData.h:215
PViewData::setFileName
virtual void setFileName(const std::string &val)
Definition: PViewData.h:75
PViewData::getAdaptiveData
adaptiveData * getAdaptiveData()
Definition: PViewData.h:238
PViewData::searchScalarClosest
bool searchScalarClosest(double x, double y, double z, double &distance, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:338
PViewData::isListBased
virtual bool isListBased()
Definition: PViewData.h:285
PViewData::getDirty
virtual bool getDirty()
Definition: PViewData.h:62
PViewData::initAdaptiveData
void initAdaptiveData(int step, int level, double tol)
Definition: PViewData.cpp:37
PViewData::_kdtree
SPoint3KDTree * _kdtree
Definition: PViewData.h:45
fullMatrix
Definition: MElement.h:27
PViewData::_fileName
std::string _fileName
Definition: PViewData.h:36
PViewData::getDimension
virtual int getDimension(int step, int ent, int ele)
Definition: PViewData.h:134
PViewData::setInterpolationMatrices
void setInterpolationMatrices(int type, const fullMatrix< double > &coefVal, const fullMatrix< double > &expVal)
Definition: PViewData.cpp:154
PViewData::getTime
virtual double getTime(int step)
Definition: PViewData.h:94
PViewData::PViewData
PViewData()
Definition: PViewData.cpp:15
PViewData::getMax
virtual double getMax(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PViewData::tagNode
virtual void tagNode(int step, int ent, int ele, int nod, int tag)
Definition: PViewData.h:148
PViewData::getType
virtual int getType(int step, int ent, int ele)
Definition: PViewData.h:183
PViewData::getNumHexahedra
virtual int getNumHexahedra(int step=-1)
Definition: PViewData.h:120
PViewData::_pc2kdtree
SPoint3CloudAdaptor< SPoint3Cloud > _pc2kdtree
Definition: PViewData.h:44
PViewData::getBoundingBox
virtual SBoundingBox3d getBoundingBox(int step=-1)=0
PViewData::getScalarValue
void getScalarValue(int step, int ent, int ele, int nod, double &val, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)
Definition: PViewData.cpp:97
SBoundingBox3d.h
PViewData::getNumScalars
virtual int getNumScalars(int step=-1)
Definition: PViewData.h:111
adaptiveData
Definition: adaptiveData.h:640
PViewData::hasTimeStep
virtual bool hasTimeStep(int step)
Definition: PViewData.h:211
PViewData::sendToServer
virtual void sendToServer(const std::string &name)
Definition: PViewDataIO.cpp:287
SPoint3CloudAdaptor< SPoint3Cloud >
PViewData::setMin
virtual void setMin(double min)=0
GModel
Definition: GModel.h:44
PViewData::fillRemoteVertexArrays
virtual int fillRemoteVertexArrays(std::string &options)
Definition: PViewData.h:282
PViewData::initAdaptiveDataLight
void initAdaptiveDataLight(int step, int level, double tol)
Definition: PViewData.cpp:47
PViewData::getNumQuadrangles
virtual int getNumQuadrangles(int step=-1)
Definition: PViewData.h:117
PViewData::removeAllInterpolationSchemes
static void removeAllInterpolationSchemes()
Definition: PViewData.cpp:209
PViewData::addMatrixToInterpolationScheme
static void addMatrixToInterpolationScheme(const std::string &name, int type, fullMatrix< double > &mat)
Definition: PViewData.cpp:219
PViewData::destroyAdaptiveData
void destroyAdaptiveData()
Definition: PViewData.cpp:79
PViewData::getNumTrihedra
virtual int getNumTrihedra(int step=-1)
Definition: PViewData.h:123
PViewData::_dirty
bool _dirty
Definition: PViewData.h:32
PViewData::~PViewData
virtual ~PViewData()
Definition: PViewData.cpp:21
PViewData::_fileNames
std::set< std::string > _fileNames
Definition: PViewData.h:37
PViewData::getModel
virtual GModel * getModel(int step)
Definition: PViewData.cpp:136
PViewData::getNumNodes
virtual int getNumNodes(int step, int ent, int ele)
Definition: PViewData.h:137
PViewData::writeMSH
virtual bool writeMSH(const std::string &fileName, double version=2.2, bool binary=false, bool saveMesh=true, bool multipleView=false, int partitionNum=-1, bool saveInterpolationMatrices=true, bool forceNodeData=false, bool forceElementData=false)
Definition: PViewDataIO.cpp:202
PViewData::deleteInterpolationMatrices
void deleteInterpolationMatrices(int type=0)
Definition: PViewData.cpp:194
PViewData::searchScalar
bool searchScalar(double x, double y, double z, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:322
PViewData::combineTime
virtual bool combineTime(nameData &nd)
Definition: PViewData.cpp:236
PViewData::setName
virtual void setName(const std::string &val)
Definition: PViewData.h:71
PViewData::_fileIndex
int _fileIndex
Definition: PViewData.h:39
MElement
Definition: MElement.h:30
SPoint3Cloud
Definition: SPoint3KDTree.h:12
PViewData::isRemote
virtual bool isRemote()
Definition: PViewData.h:281
PViewData::removeInterpolationScheme
static void removeInterpolationScheme(const std::string &name)
Definition: PViewData.cpp:199
PViewData::getElement
virtual MElement * getElement(int step, int entity, int element)
Definition: PViewData.cpp:148
element
Definition: shapeFunctions.h:12
PViewData::_interpolationSchemes
static std::map< std::string, interpolationMatrices > _interpolationSchemes
Definition: PViewData.h:53
PViewData::setDirty
virtual void setDirty(bool val)
Definition: PViewData.h:63
PViewData::getNumVectors
virtual int getNumVectors(int step=-1)
Definition: PViewData.h:112
PViewData
Definition: PViewData.h:29
PViewData::getNumTetrahedra
virtual int getNumTetrahedra(int step=-1)
Definition: PViewData.h:119
PViewData::empty
virtual bool empty()
Definition: PViewData.cpp:85
interpolationMatrices
std::map< int, std::vector< fullMatrix< double > * > > interpolationMatrices
Definition: PViewData.h:24
PViewData::getValue
virtual void getValue(int step, int ent, int ele, int nod, int comp, double &val)
Definition: PViewData.h:163
PViewData::setXYZV
virtual void setXYZV(std::vector< double > &x, std::vector< double > &y, std::vector< double > &z, std::vector< double > &v)
Definition: PViewData.h:275
PViewData::getMemoryInMb
virtual double getMemoryInMb()
Definition: PViewData.h:288
PViewData::getNumComponents
virtual int getNumComponents(int step, int ent, int ele)
Definition: PViewData.h:152
PViewData::getListPointers
virtual void getListPointers(int N[24], std::vector< double > *V[24])
Definition: PViewDataIO.cpp:282
PViewData::getFileName
virtual std::string getFileName(int step=-1)
Definition: PViewData.h:74
PViewData::getNumTensors
virtual int getNumTensors(int step=-1)
Definition: PViewData.h:113
z
const double z
Definition: GaussQuadratureQuad.cpp:56
PViewData::setMax
virtual void setMax(double max)=0
PViewData::getNumPrisms
virtual int getNumPrisms(int step=-1)
Definition: PViewData.h:121
PViewData::_adaptive
adaptiveData * _adaptive
Definition: PViewData.h:49
PViewData::importLists
virtual void importLists(int N[24], std::vector< double > *V[24])
Definition: PViewDataIO.cpp:277
PViewData::getInterpolationSchemeName
std::string getInterpolationSchemeName()
Definition: PViewData.h:258
PViewData::getSizeInterpolationScheme
static int getSizeInterpolationScheme()
Definition: PViewData.cpp:226
PViewData::getMin
virtual double getMin(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PViewData::getFileIndex
virtual int getFileIndex()
Definition: PViewData.h:86
PViewData::setFileIndex
virtual void setFileIndex(int val)
Definition: PViewData.h:87
PViewData::getNumElements
virtual int getNumElements(int step=-1, int ent=-1)
Definition: PViewData.h:131
PViewData::searchTensorClosest
bool searchTensorClosest(double x, double y, double z, double &distance, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:418
PViewData::getNumPolyhedra
virtual int getNumPolyhedra(int step=-1)
Definition: PViewData.h:124
PViewData::getName
virtual std::string getName()
Definition: PViewData.h:70
PViewData::saveAdaptedViewForVTK
void saveAdaptedViewForVTK(const std::string &guifileName, int useDefaultName, int step, int level, double tol, int npart, bool isBinary)
Definition: PViewData.cpp:61
nanoflann::KDTreeSingleIndexAdaptor
Definition: nanoflann.hpp:746
PViewData::getString2D
virtual void getString2D(int i, int step, std::string &str, double &x, double &y, double &style)
Definition: PViewData.h:190
PViewData::setBoundingBox
virtual void setBoundingBox(SBoundingBox3d &box)=0
PViewData::hasFileName
virtual bool hasFileName(const std::string &val)
Definition: PViewData.h:80
PViewData::smooth
virtual void smooth()
Definition: PViewData.cpp:231
PViewData::fromVector
virtual bool fromVector(const std::vector< std::vector< double > > &vec)
Definition: PViewDataIO.cpp:238
PViewData::searchTensor
bool searchTensor(double x, double y, double z, double *values, int step=-1, double *size=nullptr, int qn=0, double *qx=nullptr, double *qy=nullptr, double *qz=nullptr, bool grad=false, int dim=-1)
Definition: PViewData.cpp:402
SBoundingBox3d
Definition: SBoundingBox3d.h:21
PViewData::finalize
virtual bool finalize(bool computeMinMax=true, const std::string &interpolationScheme="")
Definition: PViewData.cpp:30
PViewData::setInterpolationSchemeName
void setInterpolationSchemeName(std::string name)
Definition: PViewData.h:259
PViewData::getNumTriangles
virtual int getNumTriangles(int step=-1)
Definition: PViewData.h:116
PViewData::useGaussPoints
virtual bool useGaussPoints()
Definition: PViewData.h:221
PViewData::setXY
virtual void setXY(std::vector< double > &x, std::vector< double > &y)
Definition: PViewData.h:272
PViewData::findClosestNode
double findClosestNode(double &xn, double &yn, double &zn, int step)
Definition: PViewData.cpp:248
PViewData::skipEntity
virtual bool skipEntity(int step, int ent)
Definition: PViewData.h:206
PViewData::getEntity
virtual GEntity * getEntity(int step, int entity)
Definition: PViewData.cpp:142
PViewData::getFirstNonEmptyTimeStep
virtual int getFirstNonEmptyTimeStep(int start=0)
Definition: PViewData.h:91