12#include "gadgetconfig.h"
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"
45 unsigned char SofteningClass;
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>
66 potdata_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr) : gcomm(dptr, tptr, pptr) {}
71 for(
int k = 0; k < 3; k++)
72 in->IntPos[k] = Tp->P[i].IntPos[k];
74 in->Type = Tp->P[i].getType();
76 in->SofteningClass = Tp->P[i].getSofteningClass();
81 void out2particle(potdata_out *out,
int i,
int mode)
override
84 Tp->PS[i].u.s.u.DM_Potential = out->Potential;
86 Tp->PS[i].u.s.u.DM_Potential += out->Potential;
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
94 int ninteractions = 0;
95 MyIntPosType intpos[3] = {in->IntPos[0], in->IntPos[1], in->IntPos[2]};
105 for(
int k = 0; k < numnodes; k++)
113 no = firstnode[k].
Node;
115 if(numnodes > Tree->MaxNodes)
116 Terminate(
"numnodes=%d Tree->MaxNodes=%d Tree->NumNodes=%d no=%d", numnodes, Tree->MaxNodes, Tree->NumNodes, no);
118 no = Tree->get_nodep(no)->nextnode;
121 unsigned int shmrank = Tree->TreeSharedMem_ThisTask;
128#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
129 char flag_Q2Tensor = 0;
132 if(no < Tree->MaxPart)
134 auto P = Tree->get_Pp(no, shmrank);
137 Tp->nearest_image_intpos_to_pos(P->IntPos, intpos, dxyz.
da);
145 hmax = (h_j > h_i) ? h_j : h_i;
148 no = Tree->get_nextnodep(shmrank)[no];
150 else if(no < Tree->MaxPart + Tree->MaxNodes)
154 if(no < Tree->FirstNonTopLevelNode)
162 nop = Tree->get_nodep(no, shmrank);
183 if(dist[0] < intlen && dist[1] < intlen && dist[2] < intlen)
193 Tp->nearest_image_intpos_to_pos(nop->
s.da, intpos, dxyz.
da);
199 double len = intlen * Tp->FacIntToCoord;
200 double len2 = len * len;
203 if(len2 > r2 * theta2)
234#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
242 else if(no >= Tree->ImportedNodeOffset)
244 Terminate(
"We don't expect TreePoints here");
249 Tree->tree_export_node_threads(no, target, &Thread);
251 no = Tree->Nextnode[no - Tree->MaxNodes];
260 double rinv = (r > 0) ? 1.0 / r : 0;
262 typename gtree::gfactors gfac;
264 Tree->get_gfactors_potential(gfac, r, hmax, rinv);
266 pot -= mass * gfac.fac0;
270#if(MULTIPOLE_ORDER >= 3) || (MULTIPOLE_ORDER >= 2 && defined(EXTRAPOTTERM))
273 double g1 = gfac.fac1 * rinv;
274 double g2 = gfac.fac2 * rinv * rinv;
276 double Q2dxyz2 = Q2dxyz * dxyz;
277 double Q2trace = nop->Q2Tensor.trace();
279 pot -= 0.5 * (g1 * Q2trace + g2 * Q2dxyz2);
290 return ninteractions;
294template <
typename partset>
295void fof<partset>::subfind_potential_compute(
domain<partset> *SubDomain,
int num,
int *darg)
298 potdata_comm<gravtree<partset>,
domain<partset>, partset> commpattern{SubDomain, &FoFGravTree, Tp};
304#include "../data/simparticles.h"
305template class fof<simparticles>;
307#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
308#include "../data/lcparticles.h"
309template class fof<lcparticles>;
global_data_all_processes All
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_LOCAL_PARTICLES
#define LEVEL_ALWAYS_OPEN
int32_t MySignedIntPosType
#define BITS_FOR_POSITIONS
gravtree< simparticles > gtree
unsigned char nextnode_shmrank
vector< MyIntPosType > center
unsigned char sibling_shmrank
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
unsigned char maxsofttype
unsigned char minsofttype