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 "../data/simparticles.h"
25#include "../domain/domain.h"
26#include "../gravtree/gravtree.h"
27#include "../logs/logs.h"
28#include "../logs/timer.h"
29#include "../main/simulation.h"
30#include "../mpi_utils/mpi_utils.h"
31#include "../ngbtree/ngbtree.h"
32#include "../sort/cxxsort.h"
33#include "../system/system.h"
34#include "../time_integration/timestep.h"
36inline int sph::sph_treetimestep_evaluate_particle_node_opening_criterion(pinfo &pdat,
ngbnode *nop)
54 vleft = (-vsig + (nop->
vmin[0] - SphP->
VelPred[0]));
58 vright = (vsig + (nop->
vmax[0] - SphP->
VelPred[0]));
71 if(right[0] < 0 || left[0] > 0)
76 vleft = (-vsig + (nop->
vmin[1] - SphP->
VelPred[1]));
80 vright = (vsig + (nop->
vmax[1] - SphP->
VelPred[1]));
93 if(right[1] < 0 || left[1] > 0)
98 vleft = (-vsig + (nop->
vmin[2] - SphP->
VelPred[2]));
102 vright = (vsig + (nop->
vmax[2] - SphP->
VelPred[2]));
115 if(right[2] < 0 || left[2] > 0)
121inline void sph::sph_treetimestep_check_particle_particle_interaction(pinfo &pdat,
int p,
int p_type,
unsigned char shmrank)
123#ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
124 if(skip_actual_force_computation)
144 double dist2 = dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2];
148 double dist =
sqrt(dist2);
149 double vsig = SphP_i->
Csnd + SphP_j->
Csnd -
156 dist += 2 * SphP_i->
Hsml;
159 double dt = dist / vsig;
175 double dist2 = dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2];
179 double dist =
sqrt(dist2);
180 double vsig = SphP_i->
Csnd + SphP_j->
Csnd -
187 dist += 2 * SphP_i->
Hsml;
190 double dt = dist / vsig;
200inline void sph::sph_treetimestep_open_node(pinfo &pdat,
ngbnode *nop,
int mintopleafnode,
int committed)
210 "p=%d < 0 nop->sibling=%d nop->nextnode=%d shmrank=%d nop->sibling_shmrank=%d nop->foreigntask=%d "
211 "first_nontoplevelnode=%d",
215 unsigned char next_shmrank;
222 next_shmrank = shmrank;
255 "should not happen: p=%d MaxPart=%d MaxNodes=%d ImportedNodeOffset=%d EndOfTreePoints=%d EndOfForeignNodes=%d "
260 sph_treetimestep_interact(pdat, p, type, shmrank, mintopleafnode, committed);
263 shmrank = next_shmrank;
267inline void sph::sph_treetimestep_interact(pinfo &pdat,
int no,
char no_type,
unsigned char shmrank,
int mintopleafnode,
int committed)
271 sph_treetimestep_check_particle_particle_interaction(pdat, no, no_type, shmrank);
284 int openflag = sph_treetimestep_evaluate_particle_node_opening_criterion(pdat, nop);
293 Terminate(
"this should not happen any more");
304 int min_buffer_space =
322 D->mpi_printf(
"TIMESTEP-TREEWALK: Begin\n");
392 int min_buffer_space =
394 if(min_buffer_space >= committed)
399 int mintopleaf =
WorkStack[item].MinTopLeafNode;
403 get_pinfo(target, pdat);
417 Terminate(
"item=%d: no=%d now we should be able to open it!", item, no);
420 sph_treetimestep_open_node(pdat, nop, mintopleaf, committed);
428 Terminate(
"Can't even process a single particle");
454 MPI_Allreduce(MPI_IN_PLACE, &max_ncycles, 1, MPI_INT, MPI_MAX,
D->Communicator);
466 D->mpi_printf(
"TIMESTEP-TREEWALK: took %g sec, max_ncycles = %d, part/sec = %g\n",
Logs.
timediff(ta, tb), max_ncycles,
global_data_all_processes All
void nearest_image_intpos_to_pos(const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
MySignedIntPosType pos_to_signedintpos(T posdiff)
MyIntPosType nearest_image_intpos_to_intpos_Y(const MyIntPosType a, const MyIntPosType b)
MyIntPosType nearest_image_intpos_to_intpos_Z(const MyIntPosType a, const MyIntPosType b)
MyIntPosType nearest_image_intpos_to_intpos_X(const MyIntPosType a, const MyIntPosType b)
double timediff(double t0, double t1)
int * GetNodeIDForSimulCommRank
int drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone=false)
This function drifts a particle i to time1.
TimeBinData TimeBinsHydro
void tree_based_timesteps(void)
long long sum_NumForeignPoints
workstack_data * WorkStack
sph_particle_data * get_SphPp(int n, unsigned char shmrank)
void cleanup_shared_memory_access(void)
domain< simparticles > * D
pdata get_Pp(int n, unsigned char shmrank)
static bool compare_workstack(const workstack_data &a, const workstack_data &b)
int * get_nextnodep(unsigned char shmrank)
long long sum_NumForeignNodes
void tree_fetch_foreign_nodes(enum ftype fetch_type)
void prepare_shared_memory_access(void)
foreign_sphpoint_data * Foreign_Points
fetch_data * StackToFetch
int AllocWorkStackBaseHigh
ngbnode * get_nodep(int no)
void tree_add_to_fetch_stack(ngbnode *nop, int nodetoopen, unsigned char shmrank)
void tree_add_to_work_stack(int target, int no, unsigned char shmrank, int mintopleafnode)
void tree_initialize_leaf_node_access_info(void)
int AllocWorkStackBaseLow
double mycxxsort(T *begin, T *end, Tcomp comp)
#define LEVEL_ALWAYS_OPEN
int32_t MySignedIntPosType
#define TREE_MIN_WORKSTACK_SIZE
#define NODE_TYPE_FETCHED_NODE
#define TREE_EXPECTED_CYCLES
#define NODE_TYPE_FETCHED_PARTICLE
#define NODE_TYPE_LOCAL_NODE
#define NODE_TYPE_LOCAL_PARTICLE
#define TIMER_START(counter)
Starts the timer counter.
#define TIMER_STOP(counter)
Stops the timer counter.
#define TIMER_STOPSTART(stop, start)
Stops the timer 'stop' and starts the timer 'start'.
long long GlobalNActiveParticles
unsigned char nextnode_shmrank
std::atomic< unsigned char > cannot_be_opened_locally
vector< MyIntPosType > center
unsigned char sibling_shmrank
unsigned char Nextnode_shmrank
sph_particle_data_hydrocore SphCore
double cf_atime2_hubble_a
std::atomic< integertime > Ti_Current
MySignedIntPosType center_offset_max[3]
MySignedIntPosType center_offset_min[3]
void drift_node(integertime time1, simparticles *Sp)
integertime get_Ti_Current(void)
unsigned char getType(void)
#define TREE_NUM_BEFORE_NODESPLIT