MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
treewalk.h
Go to the documentation of this file.
1 #ifndef _EVALUATOR_H_
2 #define _EVALUATOR_H_
3 
4 #include <stdint.h>
5 #include "utils/paramset.h"
6 #include "forcetree.h"
7 
8 #define NODELISTLENGTH 8
9 
13 };
14 
18 };
19 
20 typedef struct TreeWalk TreeWalk;
21 
22 typedef struct {
23  double Pos[3];
24 #ifdef DEBUG
25  MyIDType ID;
26 #endif
27  int NodeList[NODELISTLENGTH];
29 
30 typedef struct {
31 #ifdef DEBUG
32  MyIDType ID;
33 #endif
35 
36 typedef struct {
37  enum NgbTreeFindSymmetric symmetric;
38  int mask;
39  double Hsml;
40  double dist[3];
41  double r2;
42  double r;
43  int other;
45 
46 typedef struct {
48 
49  int mode; /* 0 for Primary, 1 for Secondary */
50  int target; /* defined only for primary (mode == 0) */
51 
52  /* Thread local export variables*/
53  size_t Nexport;
54  size_t BunchSize;
55  /* Number of entries in the export table for this particle*/
57  int *exportflag;
59  size_t *exportindex;
61 
62  int * ngblist;
63  int64_t Ninteractions;
64  int64_t Nnodesinlist;
65  int64_t Nlist;
67 
69 
71 
72 typedef int (*TreeWalkHasWorkFunction) (const int i, TreeWalk * tw);
73 typedef void (*TreeWalkProcessFunction) (const int i, TreeWalk * tw);
74 
75 typedef void (*TreeWalkFillQueryFunction)(const int j, TreeWalkQueryBase * query, TreeWalk * tw);
76 typedef void (*TreeWalkReduceResultFunction)(const int j, TreeWalkResultBase * result, const enum TreeWalkReduceMode mode, TreeWalk * tw);
77 
82 };
83 
84 struct TreeWalk {
85  void * priv;
86 
87  /* A pointer to the force tree structure to walk.*/
88  const ForceTree * tree;
89 
90  /* name of the evaluator (used in printing messages) */
91  const char * ev_label;
92 
93  enum TreeWalkType type;
94 
98 
99  TreeWalkVisitFunction visit; /* Function to be called between a tree node and a particle */
100  TreeWalkHasWorkFunction haswork; /* Is the particle part of this interaction? */
101  TreeWalkFillQueryFunction fill; /* Copy the useful attributes of a particle to a query */
102  TreeWalkReduceResultFunction reduce; /* Reduce a partial result to the local particle storage */
103  TreeWalkNgbIterFunction ngbiter; /* called for each pair of particles if visit is set to ngbiter */
104  TreeWalkProcessFunction postprocess; /* postprocess finalizes quantities for each particle, e.g. divide the normalization */
105  TreeWalkProcessFunction preprocess; /* Preprocess initializes quantities for each particle */
106  int NTask; /*Number of MPI tasks*/
107  int64_t NThread; /*Number of OpenMP threads*/
108 
109  /* Unlike in Gadget-3, when exporting we now always send tree branches.*/
110  char * dataget;
111  char * dataresult;
112 
113  /* The metal return alters neighbours,
114  * which cannot be evaluated twice.
115  * If repeatdisallowd is true, we allocate memory
116  * to keep track of the evaluated particles.*/
118  char * evaluated;
119 
120  /* performance metrics */
121  double timewait1;
122  double timewait2;
123  double timecomp1;
124  double timecomp2;
125  double timecomp3;
128  /* For secondary tree walks this stores the
129  * total number of pseudo-particles in all
130  * node lists of exported particles.*/
131  int64_t Nnodesinlist;
132  /* Stores the total number of node lists created for all exported particles.
133  * Used to find the average number of nodes in each nodelist.*/
134  int64_t Nlist;
135  /* Total number of exported particles
136  * (Nexport is only the exported particles in the current export buffer). */
137  int64_t Nexport_sum;
138  /* Number of times we filled up our export buffer*/
139  int64_t Nexportfull;
140  /* Number of times we needed to re-run the treewalk.
141  * Convenience variable for density. */
142  int64_t Niteration;
143 
144  /* internal flags*/
145  /* Number of particles marked for export to another processor*/
146  size_t Nexport;
147  /* Number of particles exported to this processor*/
148  size_t Nimport;
149  /* Flags that our export buffer is full*/
151  /* Number of particles we can fit into the export buffer*/
152  size_t BunchSize;
153  /* List of neighbour candidates.*/
154  int *Ngblist;
155  /* Flag not allocating nighbour list*/
157  /* Index into WorkSet to start iteration.
158  * Will be !=0 if the export buffer fills up*/
159  int64_t WorkSetStart;
160  /* The list of particles to work on. May be NULL, in which case all particles are used.*/
161  int * WorkSet;
162  /* Size of the workset list*/
163  int64_t WorkSetSize;
164  /*Did we use the active_set array as the WorkSet?*/
166  /* Redo counters and queues*/
167  size_t *NPLeft;
168  int **NPRedo;
170  /* Max and min arrays for each iteration of the count*/
171  double * maxnumngb;
172  double * minnumngb;
173 };
174 
175 /*Initialise treewalk parameters on first run*/
177 
178 /* Do the distributed tree walking. Warning: as this is a threaded treewalk,
179  * it may call tw->visit on particles more than once and in a noneterministic order.
180  * Your module should behave correctly in this case! */
181 void treewalk_run(TreeWalk * tw, int * active_set, size_t size);
182 
184  TreeWalkResultBase * O,
185  LocalTreeWalk * lv);
186 
187 /*returns -1 if the buffer is full */
188 int treewalk_export_particle(LocalTreeWalk * lv, int no);
189 #define TREEWALK_REDUCE(A, B) (A) = (mode==TREEWALK_PRIMARY)?(B):((A) + (B))
190 
191 /*****
192  * Variant of ngbiter that doesn't use the Ngblist.
193  * The ngblist is generally preferred for memory locality reasons and
194  * to avoid particles being partially evaluated
195  * twice if the buffer fills up. Use this variant if the evaluation
196  * wants to change the search radius, such as for knn algorithms
197  * or some density code. Don't use it if the treewalk modifies other particles.
198  * */
200 
201 #define MAXITER 400
202 
203 /* This function does treewalk_run in a loop, allocating a queue to allow some particles to be redone.
204  * This loop is used primarily in density estimation.*/
205 void treewalk_do_hsml_loop(TreeWalk * tw, int * queue, int64_t queuesize, int update_hsml);
206 
207 /* This function find the closest index in the multi-evaluation list of hsml and numNgb, update left and right bound, and return the new hsml */
208 double ngb_narrow_down(double *right, double *left, const double *radius, const double *numNgb, int maxcmpt, int desnumngb, int *closeidx, double BoxSize);
209 
210 void treewalk_build_queue(TreeWalk * tw, int * active_set, const size_t size, int may_have_garbage);
211 
212 #endif
size_t Nexport
Definition: treewalk.h:53
int * exportnodecount
Definition: treewalk.h:58
size_t NThisParticleExport
Definition: treewalk.h:56
size_t BunchSize
Definition: treewalk.h:54
int64_t Nlist
Definition: treewalk.h:65
int * exportflag
Definition: treewalk.h:57
int * ngblist
Definition: treewalk.h:62
int64_t Ninteractions
Definition: treewalk.h:63
size_t DataIndexOffset
Definition: treewalk.h:60
int64_t Nnodesinlist
Definition: treewalk.h:64
size_t * exportindex
Definition: treewalk.h:59
TreeWalk * tw
Definition: treewalk.h:47
int64_t Niteration
Definition: treewalk.h:142
double * minnumngb
Definition: treewalk.h:172
enum TreeWalkType type
Definition: treewalk.h:93
int64_t NThread
Definition: treewalk.h:107
char * dataresult
Definition: treewalk.h:111
void * priv
Definition: treewalk.h:85
TreeWalkHasWorkFunction haswork
Definition: treewalk.h:100
int64_t WorkSetSize
Definition: treewalk.h:163
size_t ngbiter_type_elsize
Definition: treewalk.h:97
TreeWalkProcessFunction preprocess
Definition: treewalk.h:105
double timecomp2
Definition: treewalk.h:124
double * maxnumngb
Definition: treewalk.h:171
int64_t Nlist
Definition: treewalk.h:134
double timecomp3
Definition: treewalk.h:125
int ** NPRedo
Definition: treewalk.h:168
char * evaluated
Definition: treewalk.h:118
char * dataget
Definition: treewalk.h:110
int64_t Nexport_sum
Definition: treewalk.h:137
size_t BunchSize
Definition: treewalk.h:152
double timecommsumm1
Definition: treewalk.h:126
const ForceTree * tree
Definition: treewalk.h:88
int BufferFullFlag
Definition: treewalk.h:150
TreeWalkProcessFunction postprocess
Definition: treewalk.h:104
int64_t WorkSetStart
Definition: treewalk.h:159
int NTask
Definition: treewalk.h:106
TreeWalkNgbIterFunction ngbiter
Definition: treewalk.h:103
int64_t Nexportfull
Definition: treewalk.h:139
int work_set_stolen_from_active
Definition: treewalk.h:165
TreeWalkReduceResultFunction reduce
Definition: treewalk.h:102
const char * ev_label
Definition: treewalk.h:91
double timewait2
Definition: treewalk.h:122
int * WorkSet
Definition: treewalk.h:161
int NoNgblist
Definition: treewalk.h:156
TreeWalkVisitFunction visit
Definition: treewalk.h:99
int repeatdisallowed
Definition: treewalk.h:117
size_t Redo_thread_alloc
Definition: treewalk.h:169
int64_t Nnodesinlist
Definition: treewalk.h:131
size_t result_type_elsize
Definition: treewalk.h:96
TreeWalkFillQueryFunction fill
Definition: treewalk.h:101
size_t * NPLeft
Definition: treewalk.h:167
double timecomp1
Definition: treewalk.h:123
size_t Nexport
Definition: treewalk.h:146
int * Ngblist
Definition: treewalk.h:154
double timewait1
Definition: treewalk.h:121
double timecommsumm2
Definition: treewalk.h:127
size_t Nimport
Definition: treewalk.h:148
size_t query_type_elsize
Definition: treewalk.h:95
int(* TreeWalkVisitFunction)(TreeWalkQueryBase *input, TreeWalkResultBase *output, LocalTreeWalk *lv)
Definition: treewalk.h:68
int treewalk_visit_ngbiter(TreeWalkQueryBase *I, TreeWalkResultBase *O, LocalTreeWalk *lv)
Definition: treewalk.c:924
int treewalk_visit_nolist_ngbiter(TreeWalkQueryBase *I, TreeWalkResultBase *O, LocalTreeWalk *lv)
Definition: treewalk.c:1143
void treewalk_run(TreeWalk *tw, int *active_set, size_t size)
Definition: treewalk.c:619
void(* TreeWalkReduceResultFunction)(const int j, TreeWalkResultBase *result, const enum TreeWalkReduceMode mode, TreeWalk *tw)
Definition: treewalk.h:76
TreeWalkType
Definition: treewalk.h:78
@ TREEWALK_ACTIVE
Definition: treewalk.h:79
@ TREEWALK_SPLIT
Definition: treewalk.h:81
@ TREEWALK_ALL
Definition: treewalk.h:80
#define NODELISTLENGTH
Definition: treewalk.h:8
void treewalk_build_queue(TreeWalk *tw, int *active_set, const size_t size, int may_have_garbage)
Definition: treewalk.c:394
TreeWalkReduceMode
Definition: treewalk.h:15
@ TREEWALK_PRIMARY
Definition: treewalk.h:16
@ TREEWALK_GHOSTS
Definition: treewalk.h:17
double ngb_narrow_down(double *right, double *left, const double *radius, const double *numNgb, int maxcmpt, int desnumngb, int *closeidx, double BoxSize)
Definition: treewalk.c:1366
void treewalk_do_hsml_loop(TreeWalk *tw, int *queue, int64_t queuesize, int update_hsml)
Definition: treewalk.c:1254
void(* TreeWalkProcessFunction)(const int i, TreeWalk *tw)
Definition: treewalk.h:73
void(* TreeWalkNgbIterFunction)(TreeWalkQueryBase *input, TreeWalkResultBase *output, TreeWalkNgbIterBase *iter, LocalTreeWalk *lv)
Definition: treewalk.h:70
void(* TreeWalkFillQueryFunction)(const int j, TreeWalkQueryBase *query, TreeWalk *tw)
Definition: treewalk.h:75
int treewalk_export_particle(LocalTreeWalk *lv, int no)
Definition: treewalk.c:567
void set_treewalk_params(ParameterSet *ps)
Definition: treewalk.c:55
int(* TreeWalkHasWorkFunction)(const int i, TreeWalk *tw)
Definition: treewalk.h:72
NgbTreeFindSymmetric
Definition: treewalk.h:10
@ NGB_TREEFIND_ASYMMETRIC
Definition: treewalk.h:12
@ NGB_TREEFIND_SYMMETRIC
Definition: treewalk.h:11
uint64_t MyIDType
Definition: types.h:10