GADGET-4
subfind_treepotential.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
16#include <mpi.h>
17#include <algorithm>
18#include <cmath>
19#include <cstdio>
20#include <cstdlib>
21#include <cstring>
22
23#include "../data/allvars.h"
24#include "../data/dtypes.h"
25#include "../data/mymalloc.h"
26#include "../domain/domain.h"
27#include "../fof/fof.h"
28#include "../gravtree/gravtree.h"
29#include "../logs/timer.h"
30#include "../main/simulation.h"
31#include "../mpi_utils/generic_comm.h"
32#include "../mpi_utils/mpi_utils.h"
33#include "../sort/peano.h"
34#include "../subfind/subfind.h"
35#include "../system/system.h"
36
39struct potdata_in : data_in_generic
40{
41 MyIntPosType IntPos[3];
42
43 unsigned char Type;
44#if NSOFTCLASSES > 1
45 unsigned char SofteningClass;
46#endif
47};
48
49struct potdata_out
50{
51 MyFloat Potential;
52};
53
54template <typename T_tree, typename T_domain, typename T_partset>
55class potdata_comm : public generic_comm<potdata_in, potdata_out, T_tree, T_domain, T_partset>
56{
57 public:
60 using gcomm::D;
61 using gcomm::Thread;
62 using gcomm::Tp; // This makes sure that we can access Tp from the base class without having to use "this->Tp"
63 using gcomm::Tree;
64
65 /* need to call the base class constructor explicitly */
66 potdata_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr) : gcomm(dptr, tptr, pptr) {}
67
68 /* routine that fills the relevant particle/cell data into the input structure defined above */
69 void particle2in(potdata_in *in, int i) override
70 {
71 for(int k = 0; k < 3; k++)
72 in->IntPos[k] = Tp->P[i].IntPos[k];
73
74 in->Type = Tp->P[i].getType();
75#if NSOFTCLASSES > 1
76 in->SofteningClass = Tp->P[i].getSofteningClass();
77#endif
78 }
79
80 /* routine to store or combine result data */
81 void out2particle(potdata_out *out, int i, int mode) override
82 {
83 if(mode == MODE_LOCAL_PARTICLES) /* initial store */
84 Tp->PS[i].u.s.u.DM_Potential = out->Potential;
85 else /* combine */
86 Tp->PS[i].u.s.u.DM_Potential += out->Potential;
87 }
88
89 int evaluate(int target, int mode, int thread_id, int action, potdata_in *in, int numnodes, node_info *firstnode,
90 potdata_out &out) override
91 {
92 gravnode *nop = NULL;
93 double pot = 0.0;
94 int ninteractions = 0;
95 MyIntPosType intpos[3] = {in->IntPos[0], in->IntPos[1], in->IntPos[2]};
96
97 double theta2 = All.ErrTolTheta * All.ErrTolTheta;
98
99#if NSOFTCLASSES > 1
100 double hmax, h_i = All.ForceSoftening[in->SofteningClass];
101#else
102 double hmax = All.ForceSoftening[0];
103#endif
104
105 for(int k = 0; k < numnodes; k++)
106 {
107 int no;
108
109 if(mode == 0)
110 no = Tree->MaxPart; /* root node */
111 else
112 {
113 no = firstnode[k].Node;
114
115 if(numnodes > Tree->MaxNodes)
116 Terminate("numnodes=%d Tree->MaxNodes=%d Tree->NumNodes=%d no=%d", numnodes, Tree->MaxNodes, Tree->NumNodes, no);
117
118 no = Tree->get_nodep(no)->nextnode; /* open it */
119 }
120
121 unsigned int shmrank = Tree->TreeSharedMem_ThisTask;
122
123 while(no >= 0)
124 {
125 vector<double> dxyz;
126 double r2, mass;
127
128#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
129 char flag_Q2Tensor = 0;
130#endif
131
132 if(no < Tree->MaxPart) /* single particle */
133 {
134 auto P = Tree->get_Pp(no, shmrank);
135
136 /* convert the integer distance to floating point */
137 Tp->nearest_image_intpos_to_pos(P->IntPos, intpos, dxyz.da);
138
139 r2 = dxyz.r2();
140
141 mass = P->getMass();
142
143#if NSOFTCLASSES > 1
144 double h_j = All.ForceSoftening[Tp->P[no].getSofteningClass()];
145 hmax = (h_j > h_i) ? h_j : h_i;
146#endif
147
148 no = Tree->get_nextnodep(shmrank)[no]; /* note: here shmrank cannot change */
149 }
150 else if(no < Tree->MaxPart + Tree->MaxNodes) /* we have an internal node */
151 {
152 if(mode == 1)
153 {
154 if(no < Tree->FirstNonTopLevelNode) /* we reached a top-level node again, which means that we are done with the
155 branch */
156 {
157 no = -1;
158 continue;
159 }
160 }
161
162 nop = Tree->get_nodep(no, shmrank);
163
164 if(nop->level <= LEVEL_ALWAYS_OPEN)
165 {
166 /* we always open the root node (its full node length couldn't be stored in the integer type */
167 no = nop->nextnode;
168 shmrank = nop->nextnode_shmrank;
169 continue;
170 }
171
172 MyIntPosType halflen = ((MyIntPosType)1) << ((BITS_FOR_POSITIONS - 1) - nop->level);
173 MyIntPosType intlen = halflen << 1;
174
175 {
176 /* check whether we lie very close to the cell, and if yes, open it */
177 MyIntPosType dist[3] = {nop->center[0] - intpos[0], nop->center[1] - intpos[1], nop->center[2] - intpos[2]};
178
179 dist[0] = (((MySignedIntPosType)dist[0]) >= 0) ? dist[0] : -dist[0];
180 dist[1] = (((MySignedIntPosType)dist[1]) >= 0) ? dist[1] : -dist[1];
181 dist[2] = (((MySignedIntPosType)dist[2]) >= 0) ? dist[2] : -dist[2];
182
183 if(dist[0] < intlen && dist[1] < intlen && dist[2] < intlen)
184 {
185 /* open cell */
186 no = nop->nextnode;
187 shmrank = nop->nextnode_shmrank;
188 continue;
189 }
190 }
191
192 /* converts the integer distance to floating point */
193 Tp->nearest_image_intpos_to_pos(nop->s.da, intpos, dxyz.da);
194
195 r2 = dxyz.r2();
196
197 mass = nop->mass;
198
199 double len = intlen * Tp->FacIntToCoord;
200 double len2 = len * len;
201
202 /* check Barnes-Hut opening criterion */
203 if(len2 > r2 * theta2)
204 {
205 /* open cell */
206 no = nop->nextnode;
207 shmrank = nop->nextnode_shmrank;
208 continue;
209 }
210
211#if NSOFTCLASSES > 1
212 double h_j = All.ForceSoftening[nop->maxsofttype];
213
214 if(h_j > h_i)
215 {
216 if(r2 < h_j * h_j)
217 {
219 {
220 /* open cell */
221 no = nop->nextnode;
222 shmrank = nop->nextnode_shmrank;
223 continue;
224 }
225 }
226 hmax = h_j;
227 }
228 else
229 hmax = h_i;
230#endif
231
232 /* ok, node can be used */
233
234#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
235 /* will need to account for quadrupole tensor of node */
236 flag_Q2Tensor = 1;
237#endif
238
239 no = nop->sibling;
240 shmrank = nop->sibling_shmrank;
241 }
242 else if(no >= Tree->ImportedNodeOffset) /* point from imported nodelist */
243 {
244 Terminate("We don't expect TreePoints here");
245 }
246 else /* pseudo particle */
247 {
248 if(mode == MODE_LOCAL_PARTICLES)
249 Tree->tree_export_node_threads(no, target, &Thread);
250
251 no = Tree->Nextnode[no - Tree->MaxNodes];
252 continue;
253 }
254
255 /* now evaluate the multipole moment */
256 if(mass)
257 {
258 double r = sqrt(r2);
259
260 double rinv = (r > 0) ? 1.0 / r : 0;
261
262 typename gtree::gfactors gfac;
263
264 Tree->get_gfactors_potential(gfac, r, hmax, rinv);
265
266 pot -= mass * gfac.fac0;
267
268 ninteractions++;
269
270#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
271 if(flag_Q2Tensor)
272 {
273 double g1 = gfac.fac1 * rinv;
274 double g2 = gfac.fac2 * rinv * rinv;
275 vector<MyDouble> Q2dxyz = nop->Q2Tensor * dxyz;
276 double Q2dxyz2 = Q2dxyz * dxyz;
277 double Q2trace = nop->Q2Tensor.trace();
278
279 pot -= 0.5 * (g1 * Q2trace + g2 * Q2dxyz2); // quadrupole potential
280 }
281#endif
282 }
283 }
284 }
285
286 /* now store the result */
287
288 out.Potential = pot;
289
290 return ninteractions;
291 }
292};
293
294template <typename partset>
295void fof<partset>::subfind_potential_compute(domain<partset> *SubDomain, int num, int *darg)
296{
297 /* create an object for handling the communication */
298 potdata_comm<gravtree<partset>, domain<partset>, partset> commpattern{SubDomain, &FoFGravTree, Tp};
299
300 commpattern.execute(num, darg, MODE_DEFAULT);
301}
302
303/* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
304#include "../data/simparticles.h"
305template class fof<simparticles>;
306
307#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
308#include "../data/lcparticles.h"
309template class fof<lcparticles>;
310#endif
311
312#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
T r2(void)
Definition: symtensors.h:187
T da[3]
Definition: symtensors.h:106
#define MODE_DEFAULT
Definition: constants.h:23
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
#define LEVEL_ALWAYS_OPEN
Definition: dtypes.h:383
int32_t MySignedIntPosType
Definition: dtypes.h:36
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
#define Terminate(...)
Definition: macros.h:15
expr sqrt(half arg)
Definition: half.hpp:2777
gravtree< simparticles > gtree
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
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
Definition: allvars.h:258
MyDouble mass
Definition: gravtree.h:65
vector< MyIntPosType > s
Definition: gravtree.h:66
unsigned char maxsofttype
Definition: gravtree.h:83
unsigned char minsofttype
Definition: gravtree.h:84
Definition: tree.h:77
int Node
Definition: tree.h:78