12#include "gadgetconfig.h"
20#include "../data/allvars.h"
21#include "../data/dtypes.h"
22#include "../data/intposconvert.h"
23#include "../data/mymalloc.h"
24#include "../domain/domain.h"
25#include "../gravtree/gravtree.h"
26#include "../logs/logs.h"
27#include "../logs/timer.h"
28#include "../main/simulation.h"
29#include "../mpi_utils/mpi_utils.h"
31#include "../system/system.h"
32#include "../time_integration/timestep.h"
34#ifdef ALLOW_DIRECT_SUMMATION
59 unsigned char SofteningClass;
61#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
62 unsigned char InsideOutsideFlag : 1;
65 directdata *DirectDataIn, *DirectDataAll;
74 accdata *DirectAccOut, *DirectAccIn;
84 for(
int k = 0; k < 3; k++)
85 DirectDataIn[nforces].IntPos[k] = Sp->
P[i].
IntPos[k];
87 DirectDataIn[nforces].Mass = Sp->
P[i].
getMass();
89 DirectDataIn[nforces].Type = Sp->
P[i].
getType();
93#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
94 DirectDataIn[nforces].InsideOutsideFlag = Sp->
P[i].InsideOutsideFlag;
99 MPI_Allgather(&nforces, 1, MPI_INT,
Recv_count, 1, MPI_INT,
D->Communicator);
104 for(
int j = 0; j <
D->NTask; j++)
112 DirectDataAll = (directdata *)
Mem.mymalloc(
"DirectDataAll", nimport *
sizeof(directdata));
114 for(
int j = 0; j <
D->NTask; j++)
127 DirectAccOut = (accdata *)
Mem.mymalloc(
"DirectDataOut", count *
sizeof(accdata));
130 for(
int i = 0; i < count; i++)
132 int target = i + first;
146 for(
int j = 0; j < nimport; j++)
153 double hmax = (h_j > h_i) ? h_j : h_i;
159 double r2 = dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2];
161 double mass = DirectDataAll[j].Mass;
167 double rinv = (r > 0) ? 1.0 / r : 0;
173#if defined(PLACEHIGHRESREGION)
176 if(DirectDataAll[j].InsideOutsideFlag ==
FLAG_INSIDE && DirectDataAll[target].InsideOutsideFlag ==
FLAG_INSIDE)
182 if(modify_gfactors_pm_monopole(gfac, r, rinv, mfp))
189 pot -= mass * gfac.fac0;
191 acc -= (mass * gfac.fac1 * rinv) * dxyz;
201 pot += mass * ew.
D0phi;
203 acc += mass * ew.
D1phi;
207 DirectAccOut[result_idx].Acc[0] = acc[0];
208 DirectAccOut[result_idx].Acc[1] = acc[1];
209 DirectAccOut[result_idx].Acc[2] = acc[2];
211 DirectAccOut[result_idx].Potential = pot;
217 DirectAccIn = (accdata *)
Mem.mymalloc(
"DirectDataIn", nforces *
sizeof(accdata));
219 MPI_Request *requests = (MPI_Request *)
Mem.mymalloc_movable(&requests,
"requests", 2 *
D->NTask *
sizeof(MPI_Request));
224 int send_first, send_count;
227 while(recvTask < D->NTask && sendTask < D->NTask)
229 while(send_first + send_count <
Recv_offset[recvTask])
231 if(sendTask >=
D->NTask - 1)
232 Terminate(
"sendTask >= NTask recvTask=%d sendTask=%d", recvTask, sendTask);
240 if(recvTask >=
D->NTask - 1)
241 Terminate(
"recvTask >= NTask recvTask=%d sendTask=%d", recvTask, sendTask);
246 int start = std::max<int>(
Recv_offset[recvTask], send_first);
249 if(next - start >= 1)
251 if(
D->ThisTask == sendTask)
252 MPI_Isend(DirectAccOut + start - send_first, (next - start) *
sizeof(accdata), MPI_BYTE, recvTask,
TAG_PDATA_SPH,
253 D->Communicator, &requests[n_requests++]);
255 if(
D->ThisTask == recvTask)
256 MPI_Irecv(DirectAccIn + start -
Recv_offset[recvTask], (next - start) *
sizeof(accdata), MPI_BYTE, sendTask,
TAG_PDATA_SPH,
257 D->Communicator, &requests[n_requests++]);
265 if(sendTask >=
D->NTask)
272 MPI_Waitall(n_requests, requests, MPI_STATUSES_IGNORE);
273 Mem.myfree(requests);
281 for(
int k = 0; k < 3; k++)
282 Sp->
P[i].
GravAccel[k] = DirectAccIn[nforces].Acc[k];
285 Sp->
P[i].Potential = DirectAccIn[nforces].Potential;
290 Mem.myfree(DirectAccIn);
291 Mem.myfree(DirectAccOut);
292 Mem.myfree(DirectDataAll);
293 Mem.myfree(DirectDataIn);
300 D->mpi_printf(
"GRAVDIRECT: force is done.\n");
306 double timedirect, sumt;
307 timedirect = tend - tstart;
309 MPI_Reduce(&timedirect, &sumt, 1, MPI_DOUBLE, MPI_SUM, 0,
D->Communicator);
313 fprintf(
Logs.
FdTimings,
"Nf=%9lld timebin=%d active part/task: avg=%g total-Nf=%lld\n",
316 fprintf(
Logs.
FdTimings,
" (direct) took=%g sec part/sec: %g ia/sec: %g\n", timedirect,
global_data_all_processes All
void ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag, ewald_data &fper)
void get_gfactors_monopole(gfactors &res, const T r, const T h_max, const T rinv)
void nearest_image_intpos_to_pos(const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
double getAllocatedBytesInMB(void)
TimeBinData TimeBinsGravity
#define TREE_ACTIVE_CUTTOFF_HIGHRES_PM
#define TREE_ACTIVE_CUTTOFF_BASE_PM
int myMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcount, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
#define TIMER_START(counter)
Starts the timer counter.
#define TIMER_STOP(counter)
Stops the timer counter.
long long GlobalNActiveParticles
long long TotNumDirectForces
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
unsigned char getSofteningClass(void)
unsigned char getType(void)
vector< MyFloat > GravAccel
void subdivide_evenly(long long N, int pieces, int index_bin, long long *first, int *count)
void myflush(FILE *fstream)