File: vtkPANLHaloFinder.h

package info (click to toggle)
paraview 5.1.2%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 221,108 kB
  • ctags: 236,092
  • sloc: cpp: 2,416,026; ansic: 190,891; python: 99,856; xml: 81,001; tcl: 46,915; yacc: 5,039; java: 4,413; perl: 3,108; sh: 1,974; lex: 1,926; f90: 748; asm: 471; pascal: 228; makefile: 198; objc: 83; fortran: 31
file content (262 lines) | stat: -rw-r--r-- 8,450 bytes parent folder | download | duplicates (2)
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
/*=========================================================================

 Program:   Visualization Toolkit
 Module:    vtkPANLHaloFinder.h

 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
 All rights reserved.
 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.

 This software is distributed WITHOUT ANY WARRANTY; without even
 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
 PURPOSE.  See the above copyright notice for more information.

 =========================================================================*/
#ifndef vtkPANLHaloFinder_h
#define vtkPANLHaloFinder_h

// .NAME vtkPANLHaloFinder.h -- Compute halos (clusters of particles)
//
// .SECTION Description
// Given an input a vtkUnstructuredGrid of points with arrays vx, vy, vz, and
// id, finds clumps of points (halos) using the cosmotools halo finder.  The
// first output is a vtkUnstructuredGrid similar to the input but with the data
// array fof_halo_tag appended to indicate which halo a point lies in.  The
// value will be -1 for points in no halo.  Optionally, subhalo finding can be
// turned on which will find subhalos of all halos over a certain size.  The
// subhalo id a point is in will be appended as the data array subhalo_tag.
//
// The second output is a summary of halo properties such as mass, center of mass,
// center (computed via one of several algorithms), and net velocity.  This
// vtkUnstructuredGrid has one point per halo.
//
// The third output is empty unless subhalo finding is turned on.  If subhalo
// finding is on, this output is similar to the second output except with data
// for each subhalo rather than each halo.  It contains one point per subhalo.

#include "vtkPVVTKExtensionsCosmoToolsModule.h" // For export macro
#include "vtkUnstructuredGridAlgorithm.h"

class vtkMultiProcessController;

class VTKPVVTKEXTENSIONSCOSMOTOOLS_EXPORT vtkPANLHaloFinder : public vtkUnstructuredGridAlgorithm
{
  vtkTypeMacro(vtkPANLHaloFinder, vtkUnstructuredGridAlgorithm)
public:
  static vtkPANLHaloFinder* New();
  void PrintSelf(ostream &os, vtkIndent indent);

  // Description:
  // Turns on/off the subhalo finder part of the algorithm
  // Default: Off
  vtkSetMacro(RunSubHaloFinder,bool)
  vtkGetMacro(RunSubHaloFinder,bool)
  vtkBooleanMacro(RunSubHaloFinder,bool)

  // Description:
  // Gets/Sets RL, the physical coordinate box size
  // Default: 256.0
  vtkSetMacro(RL,double)
  vtkGetMacro(RL,double)

  // Description:
  // Gets/Sets the distance conversion factor.  This is multiplied into all position
  // coordinates before the halo finder is run.
  // Default: 1.0
  vtkSetMacro(DistanceConvertFactor,double)
  vtkGetMacro(DistanceConvertFactor,double)

  // Description:
  // Gets/Sets the mass conversion factor.  This is multiplied into the particle mass
  // before the halo finder is run.
  // Default: 1.0
  vtkSetMacro(MassConvertFactor,double)
  vtkGetMacro(MassConvertFactor,double)

  // Description:
  // Gets/Sets the size of the ghost particle region around each process's particles
  // to exchange when creating ghost particles.
  // Default: 8.0
  vtkSetMacro(DeadSize,double)
  vtkGetMacro(DeadSize,double)

  // Description:
  // Gets/Sets the particle mass.  For input datasets that do not have mass information
  // the mass of each particle defaults to this value.
  // Default: 1.307087181e+09
  vtkSetMacro(ParticleMass,float)
  vtkGetMacro(ParticleMass,float)

  // Description:
  // Gets/Sets distance threshold for particles to be considered in the same
  // halo.  This is measured in grid units on a NP x NP x NP grid.
  // Default: 0.1679999998
  vtkSetMacro(BB,double)
  vtkGetMacro(BB,double)

  // Description:
  // Gets/Sets alpha factor.  This controls how aggressively small subhalos
  // are grown.  Alpha factor of 1.0 is the least aggressive
  // Default: 1.0
  vtkSetClampMacro(AlphaFactor,double,0.0,1.0)
  vtkGetMacro(AlphaFactor,double)

  // Description:
  // Gets/Sets beta factor.  This controlls how saddle points between
  // subhalos are treated.  Larger values allow identification of smaller
  // scale structures such as tails.
  // Default: 0.0
  vtkSetClampMacro(BetaFactor,double,0.0,1.0)
  vtkGetMacro(BetaFactor,double)

  // Description:
  // Gets/Sets NP
  // Default: 1024
  vtkSetMacro(NP,int)
  vtkGetMacro(NP,int)

  // Description:
  // Gets/Sets the minimum number of close neighbors for a halo candidate to
  // include a particle.
  // Default: 1
  vtkSetMacro(NMin,int)
  vtkGetMacro(NMin,int)

  // Description:
  // Gets/Sets the minimum number of particles required for a halo candidate to
  // be considered a halo and output
  // Default: 10000
  vtkSetMacro(PMin,int)
  vtkGetMacro(PMin,int)

  // Description:
  // Gets/Sets the minimum halo size to run the subhalo finder on.
  // Default: 10000
  vtkSetMacro(MinFOFSubhaloSize,long)
  vtkGetMacro(MinFOFSubhaloSize,long)

  // Description:
  // Gets/Sets the minimum size of a subhalo candidate
  // Default: 200
  vtkSetMacro(MinCandidateSize,int)
  vtkGetMacro(MinCandidateSize,int)

  // Description:
  // Gets/Sets NumSPHNeighbors
  // Default: 64
  vtkSetMacro(NumSPHNeighbors,int)
  vtkGetMacro(NumSPHNeighbors,int)

  // Description:
  // Gets/Sets the number of neighbors that are examined by the subhalo finder
  // to determine local clumps near the each particle
  // Default: 20
  vtkSetMacro(NumNeighbors,int)
  vtkGetMacro(NumNeighbors,int)

  enum CenterFindingType {
    NONE = 0,
    MOST_BOUND_PARTICLE = 1,
    MOST_CONNECTED_PARTICLE = 2,
    HIST_CENTER_FINDING = 3
  };

  // Description:
  // Gets/Sets the center finding method used by the halo finder once halos are
  // identified.
  // Default: NONE
  vtkSetMacro(CenterFindingMode,int)
  vtkGetMacro(CenterFindingMode,int)

  // Description:
  // Gets/Sets the smoothing length used by the center finders
  // Default: 0.0
  vtkSetMacro(SmoothingLength,double)
  vtkGetMacro(SmoothingLength,double)

  // Description:
  // Gets/Sets the OmegaDM parameter of the simulation.  Used by the center
  // finding algorithms.
  // Default: 0.26627
  vtkSetMacro(OmegaDM,double)
  vtkGetMacro(OmegaDM,double)

  // Description:
  // Gets/Sets the OmegaNU parameter of the simulation.  Used by the center
  // finding algorithms.
  // Default: 0.0
  vtkSetMacro(OmegaNU,double)
  vtkGetMacro(OmegaNU,double)

  // Description:
  // Gets/Sets the Deut parameter of the simulation.  Used by the center
  // finding algorithms.
  // Default: 0.02258
  vtkSetMacro(Deut,double)
  vtkGetMacro(Deut,double)

  // Description:
  // Gets/Sets the Hubble parameter of the simulation.  Used by the center
  // finding algorithms.
  // Default: 0.673
  vtkSetMacro(Hubble,double)
  vtkGetMacro(Hubble,double)

  // Description:
  // Gets/Sets the current redshift.  Used by the center finding algorithms.
  // Default: 0.0
  vtkSetMacro(RedShift,double)
  vtkGetMacro(RedShift,double)

protected:
  vtkPANLHaloFinder();
  virtual ~vtkPANLHaloFinder();

  virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
  virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
  int FillInputPortInformation(int port, vtkInformation *info);

  double RL;
  double DistanceConvertFactor;
  double MassConvertFactor;
  double DeadSize;
  float ParticleMass;
  double BB;
  double AlphaFactor;
  double BetaFactor;
  int NP;
  int NMin;
  int PMin;
  long MinFOFSubhaloSize;
  int MinCandidateSize;
  int NumSPHNeighbors;
  int NumNeighbors;

  bool RunSubHaloFinder;

  // Center finding parameters
  int CenterFindingMode;
  double SmoothingLength;
  double OmegaNU;
  double OmegaDM;
  double Deut;
  double Hubble;
  double RedShift;

  vtkMultiProcessController* Controller;

  class vtkInternals;
  vtkInternals* Internal;
private:
  vtkPANLHaloFinder(const vtkPANLHaloFinder&); // Not implemented
  void operator=(const vtkPANLHaloFinder&); // Not implemented

  void ExtractDataArrays(vtkUnstructuredGrid* input, vtkIdType offset);
  void DistributeInput();
  void CreateGhostParticles();
  void ExecuteHaloFinder(vtkUnstructuredGrid* allParticles, vtkUnstructuredGrid* fofProperties);
  void ExecuteSubHaloFinder(vtkUnstructuredGrid* allParticles, vtkUnstructuredGrid* subFofProperties);
  void FindCenters(vtkUnstructuredGrid* allParticles, vtkUnstructuredGrid* fofProperties);
};

#endif // vtkPANLHaloFinder_h