GADGET-4
subfind_findlinkngb.cc
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#include "gadgetconfig.h"
13
14#ifdef SUBFIND
15#ifndef SUBFIND_HBT
16
17#include <gsl/gsl_math.h>
18#include <math.h>
19#include <mpi.h>
20#include <stdlib.h>
21#include <string.h>
22#include <algorithm>
23#include <cstdio>
24
25#include "../data/allvars.h"
26#include "../data/dtypes.h"
27#include "../data/intposconvert.h"
28#include "../data/mymalloc.h"
29#include "../domain/domain.h"
30#include "../fof/fof.h"
31#include "../gravtree/gravtree.h"
32#include "../logs/timer.h"
33#include "../main/simulation.h"
34#include "../mpi_utils/generic_comm.h"
35#include "../mpi_utils/mpi_utils.h"
36#include "../sort/cxxsort.h"
37#include "../subfind/subfind.h"
38#include "../system/system.h"
39
40struct r2type
41{
42 MyFloat r2;
43 int index;
44};
45
46static bool subfind_ngb_compare_dist(const r2type &a, const r2type &b) { return a.r2 < b.r2; }
47
48static int *DM_NumNgb;
49static double *Dist2list;
50static MyFloat *Left, *Right;
51
52/* local data structure for collecting particle/cell data that is sent to other processors if needed */
53struct nearest_in : data_in_generic
54{
55 MyIntPosType IntPos[3];
56 MyFloat DM_Hsml;
57};
58
59struct nearest_out
60{
61 int Ngb;
62};
63
64template <typename T_tree, typename T_domain, typename T_partset>
65class nearest_comm : public generic_comm<nearest_in, nearest_out, T_tree, T_domain, T_partset>
66{
67 public:
69 using gcomm::D;
70 using gcomm::Thread;
71 using gcomm::Tp; // This makes sure that we can access Tp from the base class without having to use "this->Tp"
72 using gcomm::Tree;
73
74 /* need to call the base class constructor explicitly */
75 nearest_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr) : gcomm(dptr, tptr, pptr) {}
76
77 /* routine that fills the relevant particle/cell data into the input structure defined above */
78 void particle2in(nearest_in *in, int i)
79 {
80 in->IntPos[0] = Tp->P[i].IntPos[0];
81 in->IntPos[1] = Tp->P[i].IntPos[1];
82 in->IntPos[2] = Tp->P[i].IntPos[2];
83 in->DM_Hsml = Tp->PS[i].v.DM_Hsml;
84 }
85
86 /* routine to store or combine result data */
87 void out2particle(nearest_out *out, int i, int mode)
88 {
89 if(mode == MODE_LOCAL_PARTICLES) /* initial store */
90 DM_NumNgb[i] = out->Ngb;
91 else /* combine */
92 DM_NumNgb[i] += out->Ngb;
93 }
94
99 int evaluate(int target, int mode, int thread_id, int action, nearest_in *in, int numnodes, node_info *firstnode, nearest_out &out)
100 {
101 MyIntPosType *intpos = in->IntPos;
102 double hsml = in->DM_Hsml;
103 int numngb = 0;
104 int exported = 0;
105
106 for(int k = 0; k < numnodes; k++)
107 {
108 int no;
109 if(mode == MODE_LOCAL_PARTICLES)
110 {
111 no = Tree->MaxPart; /* root node */
112 }
113 else
114 {
115 no = firstnode[k].Node;
116 no = Tree->get_nodep(no)->nextnode; /* open it */
117 }
118
119 unsigned int shmrank = Tree->TreeSharedMem_ThisTask;
120
121 while(no >= 0)
122 {
123 if(no < Tree->MaxPart) /* single particle */
124 {
125 auto *P = Tree->get_Pp(no, shmrank);
126
127 no = Tree->get_nextnodep(shmrank)[no]; /* note: here shmrank cannot change */
128
129 double dxyz[3];
130 Tp->nearest_image_intpos_to_pos(P->IntPos, intpos, dxyz); /* converts the integer distance to floating point */
131
132 double h2 = hsml * hsml;
133
134 double r2 = dxyz[0] * dxyz[0];
135 if(r2 > h2)
136 continue;
137
138 r2 += dxyz[1] * dxyz[1];
139 if(r2 > h2)
140 continue;
141
142 r2 += dxyz[2] * dxyz[2];
143 if(r2 > h2)
144 continue;
145
146 if(numngb >= Tp->NumPart)
147 Terminate("numngb >= Tp->NumPart");
148
149 Dist2list[numngb++] = r2;
150 }
151 else if(no < Tree->MaxPart + Tree->MaxNodes) /* internal node */
152 {
153 if(mode == 1)
154 {
155 if(no < Tree->FirstNonTopLevelNode) /* we reached a top-level node again, which means that we are done with the
156 branch */
157 break;
158 }
159
160 gravnode *current = Tree->get_nodep(no, shmrank);
161
162 no = current->sibling; /* in case the node can be discarded */
163 shmrank = current->sibling_shmrank;
164
165 double dxyz[3];
166 Tp->nearest_image_intpos_to_pos(current->center.da, intpos,
167 dxyz); /* converts the integer distance to floating point */
168
169 double lenhalf = (((MyIntPosType)1) << (BITS_FOR_POSITIONS - 1 - current->level)) * Tp->FacIntToCoord;
170
171 double dist = hsml + lenhalf;
172
173 if(fabs(dxyz[0]) > dist)
174 continue;
175 if(fabs(dxyz[1]) > dist)
176 continue;
177 if(fabs(dxyz[2]) > dist)
178 continue;
179
180 /* now test against the minimal sphere enclosing everything */
181 dist += FACT1 * 2.0 * lenhalf;
182 if(dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2] > dist * dist)
183 continue;
184
185 no = current->nextnode; /* ok, we need to open the node */
186 shmrank = current->nextnode_shmrank;
187 }
188 else
189 {
190 /* pseudo particle */
191
192 if(mode == MODE_LOCAL_PARTICLES)
193 if(target >= 0) /* if no target is given, export will not occur */
194 {
195 exported = 1;
196
197 if(mode == MODE_LOCAL_PARTICLES)
198 Tree->tree_export_node_threads(no, target, &Thread);
199 }
200
201 no = Tree->Nextnode[no - Tree->MaxNodes];
202 }
203 }
204 }
205
206 if(mode == MODE_LOCAL_PARTICLES) /* local particle */
207 if(exported == 0) /* completely local */
208 if(numngb >= All.DesLinkNgb)
209 {
210 r2type *R2list = (r2type *)Mem.mymalloc("R2list", sizeof(r2type) * numngb);
211 for(int i = 0; i < numngb; i++)
212 {
213 R2list[i].r2 = Dist2list[i];
214 }
215
216 mycxxsort(R2list, R2list + numngb, subfind_ngb_compare_dist);
217
218 Tp->PS[target].v.DM_Hsml = sqrt(R2list[All.DesLinkNgb - 1].r2);
219 numngb = All.DesLinkNgb;
220
221 for(int i = 0; i < numngb; i++)
222 {
223 Dist2list[i] = R2list[i].r2;
224 }
225
226 Mem.myfree(R2list);
227 }
228
229 out.Ngb = numngb;
230
231 return 0;
232 }
233};
234
235template <typename partset>
236void fof<partset>::subfind_find_linkngb(domain<partset> *SubDomain, int num, int *list)
237{
238 subfind_collective_printf("SUBFIND: root-task=%d: Start find_linkngb. (%d particles on root-task)\n", ThisTask, num);
239
240 Dist2list = (double *)Mem.mymalloc("Dist2list", Tp->NumPart * sizeof(double));
241 Left = (MyFloat *)Mem.mymalloc("Left", sizeof(MyFloat) * Tp->NumPart);
242 Right = (MyFloat *)Mem.mymalloc("Right", sizeof(MyFloat) * Tp->NumPart);
243 DM_NumNgb = (int *)Mem.mymalloc_movable(&DM_NumNgb, "DM_NumNgb", sizeof(int) * Tp->NumPart);
244
245 int *targetlist = (int *)Mem.mymalloc("targetlist", num * sizeof(int));
246
247 for(int idx = 0; idx < num; idx++)
248 {
249 targetlist[idx] = list[idx]; /* to preserve the input list, we make a copy */
250
251 int i = list[idx];
252 Left[i] = Right[i] = 0;
253 }
254
255 nearest_comm<gravtree<partset>, domain<partset>, partset> commpattern{SubDomain, &FoFGravTree, Tp};
256
257 /* we will repeat the whole thing for those particles where we didn't find enough neighbours */
258 long long ntot;
259 int iter = 0;
260
261 do
262 {
263 double t0 = Logs.second();
264
265 commpattern.execute(num, targetlist, MODE_DEFAULT);
266
267 /* do final operations on results */
268 int npleft = 0;
269 for(int idx = 0; idx < num; idx++)
270 {
271 int i = targetlist[idx];
272
273 /* now check whether we had enough neighbours */
274
275 if(DM_NumNgb[i] != All.DesLinkNgb && ((Right[i] - Left[i]) > 1.0e-6 * Left[i] || Left[i] == 0 || Right[i] == 0))
276 {
277 /* need to redo this particle */
278 targetlist[npleft++] = i;
279
280 if(DM_NumNgb[i] < All.DesLinkNgb)
281 Left[i] = std::max<double>(Tp->PS[i].v.DM_Hsml, Left[i]);
282 else
283 {
284 if(Right[i] != 0)
285 {
286 if(Tp->PS[i].v.DM_Hsml < Right[i])
287 Right[i] = Tp->PS[i].v.DM_Hsml;
288 }
289 else
290 Right[i] = Tp->PS[i].v.DM_Hsml;
291 }
292
293 if(iter >= MAXITER - 10)
294 {
295 double pos[3];
296 Tp->intpos_to_pos(Tp->P[i].IntPos, pos);
297
298 printf("i=%d task=%d ID=%d DM_Hsml=%g Left=%g Right=%g Right-Left=%g\n pos=(%g|%g|%g)\n", i, ThisTask,
299 (int)Tp->P[i].ID.get(), Tp->PS[i].v.DM_Hsml, Left[i], Right[i], (double)(Right[i] - Left[i]), pos[0], pos[1],
300 pos[2]);
301 myflush(stdout);
302 }
303
304 if(Right[i] > 0 && Left[i] > 0)
305 Tp->PS[i].v.DM_Hsml = pow(0.5 * (pow(Left[i], 3) + pow(Right[i], 3)), 1.0 / 3);
306 else
307 {
308 if(Right[i] == 0 && Left[i] == 0)
309 Terminate("can't occur");
310
311 if(Right[i] == 0 && Left[i] > 0)
312 Tp->PS[i].v.DM_Hsml *= 1.26;
313
314 if(Right[i] > 0 && Left[i] == 0)
315 Tp->PS[i].v.DM_Hsml /= 1.26;
316 }
317 }
318 }
319
320 num = npleft;
321
322 sumup_large_ints(1, &npleft, &ntot, SubComm);
323
324 double t1 = Logs.second();
325
326 if(ntot > 0)
327 {
328 iter++;
329
330 if(iter > 0)
331 subfind_collective_printf(
332 "SUBFIND: root-task=%d: find linkngb iteration %d, need to repeat for %lld particles. (took %g sec)\n", ThisTask, iter,
333 ntot, Logs.timediff(t0, t1));
334
335 if(iter > MAXITER)
336 Terminate("failed to converge in neighbour iteration in density_findlinkngb()\n");
337 }
338 }
339 while(ntot > 0);
340
341 Mem.myfree(targetlist);
342 Mem.myfree(DM_NumNgb);
343 Mem.myfree(Right);
344 Mem.myfree(Left);
345
346 Mem.myfree(Dist2list);
347
348 subfind_collective_printf("SUBFIND: root-task=%d: Done with find_linkngb\n", ThisTask);
349}
350
351/* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
352#include "../data/simparticles.h"
353template class fof<simparticles>;
354
355#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
356#include "../data/lcparticles.h"
357template class fof<lcparticles>;
358#endif
359
360#endif
361#endif
global_data_all_processes All
Definition: main.cc:40
Definition: domain.h:31
virtual void out2particle(T_out *out, int i, int mode)=0
virtual void particle2in(T_in *in, int i)=0
virtual int evaluate(int target, int mode, int thread_id, int action, T_in *in, int numnodes, node_info *firstnode, T_out &out)=0
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
#define MAXITER
Definition: constants.h:305
#define MODE_DEFAULT
Definition: constants.h:23
#define FACT1
Definition: constants.h:435
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
double mycxxsort(T *begin, T *end, Tcomp comp)
Definition: cxxsort.h:39
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:15
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr pow(half base, half exp)
Definition: half.hpp:2803
expr sqrt(half arg)
Definition: half.hpp:2777
unsigned char level
Definition: tree.h:64
unsigned char nextnode_shmrank
Definition: tree.h:66
int nextnode
Definition: tree.h:57
vector< MyIntPosType > center
Definition: tree.h:51
unsigned char sibling_shmrank
Definition: tree.h:65
int sibling
Definition: tree.h:53
Definition: tree.h:77
int Node
Definition: tree.h:78
void myflush(FILE *fstream)
Definition: system.cc:125