This file is indexed.

/usr/include/gmsh/PViewData.h is in libgmsh-dev 2.8.5+dfsg-1.1+b1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
// Gmsh - Copyright (C) 1997-2014 C. Geuzaine, J.-F. Remacle
//
// See the LICENSE.txt file for license information. Please report all
// bugs and problems to the public mailing list <gmsh@geuz.org>.

#ifndef _PVIEW_DATA_H_
#define _PVIEW_DATA_H_

#include <string>
#include <vector>
#include <map>
#include <set>
#include "SBoundingBox3d.h"
#include "fullMatrix.h"

#define VAL_INF 1.e200

class adaptiveData;
class GModel;
class GEntity;
class MElement;
class nameData;
class OctreePost;

typedef std::map<int, std::vector<fullMatrix<double>*> > interpolationMatrices;

// The abstract interface to post-processing view data.
class PViewData {
 private:
  // flag to mark that the data is 'dirty' and should not be displayed
  bool _dirty;
  // name of the view
  std::string _name;
  // name of the file the data was loaded from
  std::string _fileName;
  std::set<std::string> _fileNames;
  // index of the view in the file
  int _fileIndex;
  // octree for rapid search
  OctreePost *_octree;

 protected:
  // adaptive visualization data
  adaptiveData *_adaptive;
  // interpolation matrices, indexed by the type of element
  interpolationMatrices _interpolation;
  // global map of "named" interpolation matrices
  static std::map<std::string, interpolationMatrices> _interpolationSchemes;

 public:
  PViewData();
  virtual ~PViewData();

  // get/set the dirty ("not ready for display") flag
  virtual bool getDirty(){ return _dirty; }
  virtual void setDirty(bool val){ _dirty = val; }

  // finalize the view data (compute min/max, etc.)
  virtual bool finalize(bool computeMinMax=true,
                        const std::string &interpolationScheme="");

  // get/set name
  virtual std::string getName(){ return _name; }
  virtual void setName(const std::string &val){ _name = val; }

  // get/set (the main) filename containing the data
  virtual std::string getFileName(int step=-1){ return _fileName; }
  virtual void setFileName(const std::string &val)
  {
    _fileName = val;
    _fileNames.insert(val);
  }
  virtual bool hasFileName(const std::string &val)
  {
    return _fileNames.find(val) != _fileNames.end();
  }

  // get/set index of view data in file
  virtual int getFileIndex(){ return _fileIndex; }
  virtual void setFileIndex(int val){ _fileIndex = val; }

  // get number of time steps in the data
  virtual int getNumTimeSteps() = 0;
  virtual int getFirstNonEmptyTimeStep(int start=0){ return start; }

  // get the time value associated with the step-th time step
  virtual double getTime(int step){ return 0.; }

  // get/set min/max for given step (global over all steps if step=-1)
  virtual double getMin(int step=-1, bool onlyVisible=false,
                        int forceNumComponents=0, int componentMap[9]=0) = 0;
  virtual double getMax(int step=-1, bool onlyVisible=false,
                        int forceNumComponents=0, int componentMap[9]=0) = 0;
  virtual void setMin(double min) = 0;
  virtual void setMax(double max) = 0;

  // get/set the bounding box
  virtual SBoundingBox3d getBoundingBox(int step=-1) = 0;
  virtual void setBoundingBox(SBoundingBox3d& box) = 0;

  // get the number of elements of a given type, for a given step
  virtual int getNumScalars(int step=-1){ return 0; }
  virtual int getNumVectors(int step=-1){ return 0; }
  virtual int getNumTensors(int step=-1){ return 0; }
  virtual int getNumPoints(int step=-1){ return 0; }
  virtual int getNumLines(int step=-1){ return 0; }
  virtual int getNumTriangles(int step=-1){ return 0; }
  virtual int getNumQuadrangles(int step=-1){ return 0; }
  virtual int getNumPolygons(int step=-1){ return 0; }
  virtual int getNumTetrahedra(int step=-1){ return 0; }
  virtual int getNumHexahedra(int step=-1){ return 0; }
  virtual int getNumPrisms(int step=-1){ return 0; }
  virtual int getNumPyramids(int step=-1){ return 0; }
  virtual int getNumPolyhedra(int step=-1){ return 0; }

  // return the number of geometrical entities in the view
  virtual int getNumEntities(int step=-1){ return 0; }

  // return the number of elements in the ent-th entity, or the total number of
  // elements if ent < 0
  virtual int getNumElements(int step=-1, int ent=-1){ return 0; }

  // return the geometrical dimension of the ele-th element in the ent-th entity
  virtual int getDimension(int step, int ent, int ele){ return 0; }

  // return the number of nodes of the ele-th element in the ent-th entity
  virtual int getNumNodes(int step, int ent, int ele){ return 0; }

  // get/set the coordinates and tag of the nod-th node from the ele-th element
  // in the ent-th entity (if the node has a tag, getNode returns it)
  virtual int getNode(int step, int ent, int ele, int nod,
                      double &x, double &y, double &z){ return 0; }
  virtual void setNode(int step, int ent, int ele, int nod,
                       double x, double y, double z);
  virtual void tagNode(int step, int ent, int ele, int nod, int tag){}

  // return the number of components available for the ele-th element in the
  // ent-th entity
  virtual int getNumComponents(int step, int ent, int ele){ return 0; }

  // return the number of values available for the ele-th element in the ent-th
  // entity
  virtual int getNumValues(int step, int ent, int ele){ return 0; }

  // get the idx'th value for the ele-th element in the ent-th entity
  virtual void getValue(int step, int ent, int ele, int idx, double &val){}

  // gets/set the comp-th component (at the step-th time step) associated with
  // the node-th node from the ele-th element in the ent-th entity
  virtual void getValue(int step, int ent, int ele, int nod, int comp, double &val){}
  virtual void setValue(int step, int ent, int ele, int nod, int comp, double val);

  // return a scalar value (same as value for scalars, norm for vectors, etc.)
  // associated with the node-th node from the ele-th element in the ent-th
  // entity
  void getScalarValue(int step, int ent, int ele, int nod, double &val,
                      int forceNumComponents=0, int componentMap[9]=0);

  // return the number of edges of the ele-th element in the ent-th
  // entity
  virtual int getNumEdges(int step, int ent, int ele){ return 0; }

  // return the type of the ele-th element in the ent-th entity
  virtual int getType(int step, int ent, int ele){ return 0; }

  // return the number of 2D/3D strings in the view
  virtual int getNumStrings2D(){ return 0; }
  virtual int getNumStrings3D(){ return 0; }

  // return the i-th 2D/3D string in the view
  virtual void getString2D(int i, int step, std::string &str,
                           double &x, double &y, double &style){}
  virtual void getString3D(int i, int step, std::string &str,
                           double &x, double &y, double &z, double &style){}

  // change the orientation of the ele-th element
  virtual void reverseElement(int step, int ent, int ele){}

  // check if the view is empty
  virtual bool empty();

  // check if we should skip the given entity/element
  virtual bool skipEntity(int step, int ent){ return false; }
  virtual bool skipElement(int step, int ent, int ele, bool checkVisibility=false,
                           int samplingRate=1);

  // check if the data has the given step/partition/etc.
  virtual bool hasTimeStep(int step){ return step >= 0 && step < getNumTimeSteps(); }
  virtual bool hasPartition(int step, int part){ return false; }
  virtual bool hasMultipleMeshes(){ return false; }
  virtual bool hasModel(GModel *model, int step=-1){ return false; }
  virtual bool isNodeData(){ return false; }

  // true if data is given at Gauss points (instead of vertices)
  virtual bool useGaussPoints(){ return false; }

  // initialize/destroy adaptive data
  void initAdaptiveData(int step, int level, double tol);
  void destroyAdaptiveData();

  // return the adaptive data
  adaptiveData *getAdaptiveData(){ return _adaptive; }

  // set/get the interpolation matrices for elements with type "type"
  // FIXME : too much overload :(
  void setInterpolationMatrices(int type,
                                const fullMatrix<double> &coefVal,
                                const fullMatrix<double> &expVal);
  void setInterpolationMatrices(int type,
                                const fullMatrix<double> &coefVal,
                                const fullMatrix<double> &expVal,
                                const fullMatrix<double> &coefGeo,
                                const fullMatrix<double> &expGeo);
  int getInterpolationMatrices(int type, std::vector<fullMatrix<double>*> &p);
  bool haveInterpolationMatrices(int type=0);

  // access to global interpolation schemes
  static void removeInterpolationScheme(const std::string &name);
  static void addMatrixToInterpolationScheme(const std::string &name, int type,
                                             fullMatrix<double> &mat);

  // smooth the data in the view (makes it C0)
  virtual void smooth();

  // combine time steps or elements from multiple datasets
  virtual bool combineTime(nameData &nd);
  virtual bool combineSpace(nameData &nd);

  // set simple X-Y data
  virtual void setXY(std::vector<double> &x, std::vector<double> &y){}

  // ask to fill vertex arrays remotely
  virtual bool isRemote(){ return false; }
  virtual int fillRemoteVertexArrays(std::string &options){ return 0; }

  // is the view a list-based dataset
  virtual bool isListBased(){ return false; }

  // get (approx) memry used by data in Mb
  virtual double getMemoryInMb(){ return 0; }

  // get GModel (if view supports it)
  virtual GModel *getModel(int step);

  // get GEntity (if view supports it)
  virtual GEntity *getEntity(int step, int entity);

  // get MElement (if view supports it)
  virtual MElement *getElement(int step, int entity, int element);

  // search for the value of the View at point x, y, z. Values are interpolated
  // using standard first order shape functions in the post element. If several
  // time steps are present, they are all interpolated unless time step is set
  // to a different value than -1.
  bool searchScalar(double x, double y, double z, double *values,
                    int step=-1, double *size=0, int qn=0,
                    double *qx=0, double *qy=0, double *qz=0);
  bool searchScalarWithTol(double x, double y, double z, double *values,
                           int step=-1, double *size=0, double tol=1.e-2, int qn=0,
                           double *qx=0, double *qy=0, double *qz=0);
  bool searchVector(double x, double y, double z, double *values,
                    int step=-1, double *size=0, int qn=0,
                    double *qx=0, double *qy=0, double *qz=0);
  bool searchVectorWithTol(double x, double y, double z, double *values,
                           int step=-1, double *size=0, double tol=1.e-2, int qn=0,
                           double *qx=0, double *qy=0, double *qz=0);
  bool searchTensor(double x, double y, double z, double *values,
                    int step=-1, double *size=0, int qn=0,
                    double *qx=0, double *qy=0, double *qz=0);
  bool searchTensorWithTol(double x, double y, double z, double *values,
                           int step=-1, double *size=0, double tol=1.e-2, int qn=0,
                           double *qx=0, double *qy=0, double *qz=0);

  // I/O routines
  virtual bool writeSTL(const std::string &fileName);
  virtual bool writeTXT(const std::string &fileName);
  virtual bool writePOS(const std::string &fileName, bool binary=false,
                        bool parsed=true, bool append=false);
  virtual bool writeMSH(const std::string &fileName, double version=2.2, bool binary=false,
                        bool saveMesh=true, bool multipleView=false,
                        int partitionNum=0, bool saveInterpolationMatrices=true,
                        bool forceNodeData=false);
  virtual bool writeMED(const std::string &fileName);
  virtual bool toVector(std::vector<std::vector<double> > &vec);
  virtual bool fromVector(const std::vector<std::vector<double> > &vec);
  virtual void importLists(int N[24], std::vector<double> *V[24]);
  virtual void getListPointers(int N[24], std::vector<double> *V[24]);
};

class nameData{
 public:
  std::string name;
  std::vector<int> indices;
  std::vector<PViewData*> data;
};

#endif