GADGET-4
subfind_nearesttwo.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 <math.h>
18#include <mpi.h>
19#include <stdio.h>
20#include <stdlib.h>
21#include <string.h>
22
23#include "../data/allvars.h"
24#include "../data/dtypes.h"
25#include "../data/intposconvert.h"
26#include "../data/mymalloc.h"
27#include "../domain/domain.h"
28#include "../fof/fof.h"
29#include "../gravtree/gravtree.h"
30#include "../logs/timer.h"
31#include "../main/simulation.h"
32#include "../mpi_utils/generic_comm.h"
33#include "../mpi_utils/mpi_utils.h"
34#include "../subfind/subfind.h"
35#include "../system/system.h"
36
40/* local data structure for collecting particle/cell data that is sent to other processors if needed */
41struct sngb_in : data_in_generic
42{
43 MyIntPosType IntPos[3];
44 MyIDType ID;
45 MyFloat Hsml;
46 MyFloat Density;
47 MyFloat Dist[2];
48 int Count;
49 location Index[2];
50};
51
52struct sngb_out
53{
54 MyFloat Dist[2];
55 location Index[2];
56 int Count;
57};
58
59template <typename T_tree, typename T_domain, typename T_partset>
60class sngb_comm : public generic_comm<sngb_in, sngb_out, T_tree, T_domain, T_partset>
61{
62 public:
64 using gcomm::D;
65 using gcomm::Thread;
66 using gcomm::Tp; // This makes sure that we can access Tp from the base class without having to use "this->Tp"
67 using gcomm::Tree;
68
69 /* need to call the base class constructor explicitly */
70 sngb_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr) : gcomm(dptr, tptr, pptr) {}
71
72 /* routine that fills the relevant particle/cell data into the input structure defined above */
73 void particle2in(sngb_in *in, int i) override
74 {
75 in->IntPos[0] = Tp->P[i].IntPos[0];
76 in->IntPos[1] = Tp->P[i].IntPos[1];
77 in->IntPos[2] = Tp->P[i].IntPos[2];
78
79 in->Hsml = Tp->PS[i].v.DM_Hsml;
80 in->ID = Tp->P[i].ID.get();
81 in->Density = Tp->PS[i].u.s.u.DM_Density;
82 in->Count = Tp->PS[i].nearest.count;
83 for(int k = 0; k < Tp->PS[i].nearest.count; k++)
84 {
85 in->Dist[k] = Tp->R2Loc[i].dist[k];
86 in->Index[k] = Tp->PS[i].nearest.index[k];
87 }
88 }
89
90 /* routine to store or combine result data */
91 void out2particle(sngb_out *out, int i, int mode) override
92 {
93 if(mode == MODE_LOCAL_PARTICLES) /* initial store */
94 {
95 Tp->PS[i].nearest.count = out->Count;
96
97 for(int k = 0; k < out->Count; k++)
98 {
99 Tp->R2Loc[i].dist[k] = out->Dist[k];
100 Tp->PS[i].nearest.index[k] = out->Index[k];
101 }
102 }
103 else /* combine */
104 {
105 for(int k = 0; k < out->Count; k++)
106 {
107 if(Tp->PS[i].nearest.count >= 1)
108 if(Tp->PS[i].nearest.index[0] == out->Index[k])
109 continue;
110
111 if(Tp->PS[i].nearest.count == 2)
112 if(Tp->PS[i].nearest.index[1] == out->Index[k])
113 continue;
114
115 int l;
116
117 if(Tp->PS[i].nearest.count < 2)
118 {
119 l = Tp->PS[i].nearest.count;
120 Tp->PS[i].nearest.count++;
121 }
122 else
123 {
124 l = (Tp->R2Loc[i].dist[0] > Tp->R2Loc[i].dist[1]) ? 0 : 1;
125
126 if(out->Dist[k] >= Tp->R2Loc[i].dist[l])
127 continue;
128 }
129
130 Tp->R2Loc[i].dist[l] = out->Dist[k];
131 Tp->PS[i].nearest.index[l] = out->Index[k];
132
133 if(Tp->PS[i].nearest.count == 2)
134 if(Tp->PS[i].nearest.index[0] == Tp->PS[i].nearest.index[1])
135 Terminate("this is not supposed to happen");
136 }
137 }
138 }
139
143 int evaluate(int target, int mode, int thread_id, int action, sngb_in *in, int numnodes, node_info *firstnode,
144 sngb_out &out) override
145 {
146 MyIntPosType *intpos = in->IntPos;
147 MyIDType ID = in->ID;
148 double density = in->Density;
149 double hsml = in->Hsml;
150 int count = in->Count;
151
152 location index[2];
153 double dist[2];
154 for(int k = 0; k < count; k++)
155 {
156 dist[k] = in->Dist[k];
157 index[k] = in->Index[k];
158 }
159
160 if(count == 2)
161 if(index[0] == index[1])
162 {
163 Terminate("target=%d mode=%d\n", target, mode);
164 }
165
166 count = 0;
167
168 hsml *= 1.00001; /* prevents that the most distant neighbour on the edge of the search region may not be found.
169 * (needed for consistency with serial algorithm)
170 */
171
172 for(int k = 0; k < numnodes; k++)
173 {
174 int no;
175
176 if(mode == MODE_LOCAL_PARTICLES)
177 {
178 no = Tree->MaxPart; /* root node */
179 }
180 else
181 {
182 no = firstnode[k].Node;
183 no = Tree->get_nodep(no)->nextnode; /* open it */
184 }
185
186 int shmrank = Tree->TreeSharedMem_ThisTask;
187
188 while(no >= 0)
189 {
190 if(no < Tree->MaxPart) /* single particle */
191 {
192 auto *P = Tree->get_Pp(no, shmrank);
193 subfind_data *PS = Tree->get_PSp(no, shmrank);
194
195 no = Tree->get_nextnodep(shmrank)[no]; /* note: here shmrank cannot change */
196
197 if(P->ID.get() != ID) /* exclude the target particle itself */
198 {
199 if(PS->u.s.u.DM_Density > density) /* we only need to look at neighbors that are denser */
200 {
201 /* converts the integer distance to floating point */
202 double dxyz[3];
203 Tp->nearest_image_intpos_to_pos(P->IntPos, intpos, dxyz);
204
205 double h2 = hsml * hsml;
206
207 double r2 = dxyz[0] * dxyz[0];
208 if(r2 > h2)
209 continue;
210
211 r2 += dxyz[1] * dxyz[1];
212 if(r2 > h2)
213 continue;
214
215 r2 += dxyz[2] * dxyz[2];
216 if(r2 > h2)
217 continue;
218
219 // ok, we found a particle. Only store up to two closest ones.
220
221 int task = D->ThisTask + (shmrank - Tree->TreeSharedMem_ThisTask);
222
223 if(task < 0 || task >= D->NTask)
224 Terminate("illegal task=%d D->NTask=%d", task, D->NTask);
225
226 if(count < 2)
227 {
228 dist[count] = r2;
229 index[count] = {task, PS->InvIndex}; /* note: ThisTask refers here to the Subdomain */
230 count++;
231 }
232 else
233 {
234 int k = (dist[0] > dist[1]) ? 0 : 1;
235
236 if(r2 < dist[k])
237 {
238 dist[k] = r2;
239 index[k] = {task, PS->InvIndex}; /* note: ThisTask refers here to the Subdomain */
240 }
241 }
242 }
243 }
244 }
245 else if(no < Tree->MaxPart + Tree->MaxNodes) /* internal node */
246 {
247 if(mode == 1)
248 {
249 if(no < Tree->FirstNonTopLevelNode) /* we reached a top-level node again, which means that we are done with the
250 branch */
251 {
252 break;
253 }
254 }
255
256 gravnode *current = Tree->get_nodep(no, shmrank);
257
258 no = current->sibling; /* in case the node can be discarded */
259 shmrank = current->sibling_shmrank;
260
261 double dxyz[3];
262 Tp->nearest_image_intpos_to_pos(current->center.da, intpos,
263 dxyz); /* converts the integer distance to floating point */
264
265 double lenhalf = (((MyIntPosType)1) << (BITS_FOR_POSITIONS - 1 - current->level)) * Tp->FacIntToCoord;
266
267 double dist = hsml + lenhalf;
268
269 if(fabs(dxyz[0]) > dist)
270 continue;
271 if(fabs(dxyz[1]) > dist)
272 continue;
273 if(fabs(dxyz[2]) > dist)
274 continue;
275
276 /* now test against the minimal sphere enclosing everything */
277 dist += FACT1 * 2.0 * lenhalf;
278 if(dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2] > dist * dist)
279 continue;
280
281 no = current->nextnode; /* ok, we need to open the node */
282 shmrank = current->nextnode_shmrank;
283 }
284 else if(no >= Tree->ImportedNodeOffset) /* point from imported nodelist */
285 {
286 Terminate("do not expect imported points here");
287 }
288 else /* pseudo particle */
289 {
290 if(mode == MODE_LOCAL_PARTICLES)
291 if(target >= 0) /* note: if no target is given, export will not occur */
292 Tree->tree_export_node_threads(no, target, &Thread);
293
294 no = Tree->Nextnode[no - Tree->MaxNodes];
295 }
296 }
297 }
298
299 out.Count = count;
300
301 for(int k = 0; k < count; k++)
302 {
303 out.Dist[k] = dist[k];
304 out.Index[k] = index[k];
305 }
306
307 return 0;
308 }
309};
310
311template <typename partset>
312void fof<partset>::subfind_find_nearesttwo(domain<partset> *SubDomain, int num, int *list)
313{
314 subfind_collective_printf("SUBFIND: root-task=%d: Start finding nearest two.\n", ThisTask);
315
316 for(int i = 0; i < num; i++)
317 Tp->PS[list[i]].nearest.count = 0;
318
319 /* create an object for handling the communication */
320 sngb_comm<gravtree<partset>, domain<partset>, partset> commpattern{SubDomain, &FoFGravTree, Tp};
321
322 commpattern.execute(num, list, MODE_DEFAULT);
323
324 subfind_collective_printf("SUBFIND: root-task=%d: Done with nearest two.\n", ThisTask);
325}
326
327/* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
328#include "../data/simparticles.h"
329template class fof<simparticles>;
330
331#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
332#include "../data/lcparticles.h"
333template class fof<lcparticles>;
334#endif
335
336#endif
337#endif
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
#define MODE_DEFAULT
Definition: constants.h:23
#define FACT1
Definition: constants.h:435
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
unsigned int MyIDType
Definition: dtypes.h:68
#define Terminate(...)
Definition: macros.h:15
half fabs(half arg)
Definition: half.hpp:2627
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
union subfind_data::@0 u