17#include "../data/allvars.h"
18#include "../data/constants.h"
19#include "../data/dtypes.h"
20#include "../data/intposconvert.h"
21#include "../data/macros.h"
22#include "../data/mymalloc.h"
23#include "../data/particle_data.h"
24#include "../data/sph_particle_data.h"
25#include "../main/main.h"
26#include "../mpi_utils/mpi_utils.h"
27#include "../mpi_utils/setcomm.h"
28#include "../system/system.h"
29#include "../time_integration/timestep.h"
69 memcpy(
static_cast<void *
>(Ptarget),
static_cast<void *
>(Psource),
sizeof(
particle_data));
74#if defined(LIGHTCONE_PARTICLES_GROUPS) && defined(FOF)
75 double *DistanceOrigin;
78#ifdef SUBFIND_ORPHAN_TREATMENT
90 int *Head, *Next, *Tail, *MinIDTask;
96 unsigned char Nonlocal : 2, MinIDChanged : 2, Marked : 2;
101 inline void link_two_particles(
int target,
int j)
103 if(Head[target] != Head[j])
106 if(Len[Head[target]] > Len[Head[j]])
116 Next[Tail[Head[p]]] = Head[s];
118 Tail[Head[p]] = Tail[Head[s]];
120 Len[Head[p]] += Len[Head[s]];
122 if(MinID[Head[s]].get() < MinID[Head[p]].get())
124 MinID[Head[p]] = MinID[Head[s]];
125 MinIDTask[Head[p]] = MinIDTask[Head[s]];
131 while((ss = Next[ss]) >= 0);
136 struct nearest_r2_data
145 double Asmth[2], Rcut[2];
148#if defined(PMGRID) && (!defined(PERIODIC) || defined(PLACEHIGHRESREGION))
149 double TotalMeshSize[2];
158#if defined(RANDOMIZE_DOMAINCENTER_TYPES) || defined(PLACEHIGHRESREGION)
163#ifdef PLACEHIGHRESREGION
167 center[1] - halflen - ReferenceIntPos[
HIGH_MESH][1],
168 center[2] - halflen - ReferenceIntPos[
HIGH_MESH][2]};
171 center[1] + halflen - ReferenceIntPos[
HIGH_MESH][1],
172 center[2] + halflen - ReferenceIntPos[
HIGH_MESH][2]};
187 inline int check_high_res_point_location(
MyIntPosType *intpos)
190 intpos[2] - ReferenceIntPos[
HIGH_MESH][2]};
216#ifdef REARRANGE_OPTION
219 if(a.TreeID < b.TreeID)
222 if(a.TreeID > b.TreeID)
255#ifdef PRESSURE_ENTROPY_SPH
264 for(
int i = 0; i <
NumPart; i++)
294 mpi_printf(
"ALLOCATE: Changing to MaxPart = %d\n", maxpartNew);
306 mpi_printf(
"ALLOCATE: Changing to MaxPartSph = %d\n", maxpartsphNew);
325 sprintf(buffer,
"particles_%d.dat",
ThisTask);
326 if((fd = fopen(buffer,
"w")))
328 fwrite(&
NumPart, 1,
sizeof(
int), fd);
329 for(
int i = 0; i <
NumPart; i++)
331 for(
int i = 0; i <
NumPart; i++)
332 fwrite(&
P[i].Vel[0], 3,
sizeof(
MyFloat), fd);
333 for(
int i = 0; i <
NumPart; i++)
348 printf(
"Task=%d, ID=%llu, Type=%d, TimeBinGrav=%d, TimeBinHydro=%d, Mass=%g, pos=%g|%g|%g, vel=%g|%g|%g, OldAcc=%g\n",
ThisTask,
351#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
352 printf(
"GravAccel=%g|%g|%g, GravPM=%g|%g|%g, Soft=%g, SoftClass=%d\n",
P[i].GravAccel[0],
P[i].GravAccel[1],
P[i].GravAccel[2],
356 printf(
"GravAccel=%g|%g|%g, Soft=%g, SoftType=%d\n",
P[i].GravAccel[0],
P[i].GravAccel[1],
P[i].GravAccel[2],
361 if(
P[i].getType() == 0)
363 printf(
"rho=%g, hsml=%g, entr=%g, csnd=%g\n",
SphP[i].Density,
SphP[i].Hsml,
SphP[i].Entropy,
SphP[i].get_sound_speed());
364 printf(
"ID=%llu SphP[p].CurrentMaxTiStep=%g\n", (
unsigned long long)
P[i].ID.get(),
SphP[i].
CurrentMaxTiStep);
376 for(
int i = 0; i <
NumPart; i++)
377 if(
P[i].ID.get() == ID)
384#ifdef HIERARCHICAL_GRAVITY
400#if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
404#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
405 void find_long_range_step_constraint(
void);
416#ifdef ADAPTIVE_HYDRO_SOFTENING
417 int get_softeningtype_for_hydro_particle(
int i)
432#ifdef INDIVIDUAL_GRAVITY_SOFTENING
435#if(INDIVIDUAL_GRAVITY_SOFTENING) & 2
436#error "INDIVIDUAL_GRAVITY_SOFTENING may not include particle type 1 which is used as a reference point"
439#if((INDIVIDUAL_GRAVITY_SOFTENING)&1) && defined(ADAPTIVE_HYDRO_SOFTENING)
440#error "INDIVIDUAL_GRAVITY_SOFTENING may not include particle type 0 when ADAPTIVE_HYDRO_SOFTENING is used"
443 int get_softening_type_from_mass(
double mass)
446 double eps = get_desired_softening_from_mass(mass);
474 double get_desired_softening_from_mass(
double mass)
483 void init_individual_softenings(
void)
489 for(
int i = 0; i <
NumPart; i++)
490 if(
P[i].getType() == 1)
495 if(massmin >
P[i].getMass())
498 if(massmax <
P[i].getMass())
503 MPI_Allreduce(&mass, &masstot, 1, MPI_DOUBLE, MPI_SUM,
Communicator);
505 MPI_Allreduce(MPI_IN_PLACE, &massmin, 1, MPI_DOUBLE, MPI_MIN,
Communicator);
506 MPI_Allreduce(MPI_IN_PLACE, &massmax, 1, MPI_DOUBLE, MPI_MAX,
Communicator);
508 All.AvgType1Mass = masstot / ndmtot;
510 mpi_printf(
"INIT: AvgType1Mass = %g (min=%g max=%g) Ndm1tot=%lld\n",
All.AvgType1Mass, massmin, massmax, ndmtot);
512 if(massmax > 1.00001 * massmin)
513 Terminate(
"Strange: Should use constant mass type-1 particles if INDIVIDUAL_GRAVITY_SOFTENING is used\n");
519 mpi_printf(
"INIT: For this AvgType1Mass, the mean particle spacing is %g and the assigned softening is %g\n",
523 for(
int i = 0; i <
NumPart; i++)
524 if(((1 <<
P[i].getType()) & (INDIVIDUAL_GRAVITY_SOFTENING)))
global_data_all_processes All
void intpos_to_pos(MyIntPosType *intpos, T *posdiff)
void mpi_printf(const char *fmt,...)
int get_timestep_bin(integertime ti_step)
void copy_particle(particle_data *Ptarget, particle_data *Psource)
unsigned short int MarkerValue
int get_active_index(int idx)
MyFloat get_DtHsml(int i)
void timebin_cleanup_list_of_active_particles(void)
void fill_active_gravity_list_with_all_particles(void)
integertime get_timestep_grav(int p)
void assign_hydro_timesteps(void)
void print_particle_info_from_ID(MyIDType ID)
Print information relative to a particle / cell to standard output given its ID. *.
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.
void reallocate_memory_maxpartsph(int maxpartsphNew)
integertime get_timestep_hydro(int p)
void reconstruct_timebins(void)
double get_utherm_from_entropy(int i)
MyFloat get_OldAcc(int i)
void reallocate_memory_maxpart(int maxpartNew)
void allocate_memory(void)
simparticles(MPI_Comm comm)
int test_if_grav_timestep_is_too_large(int p, int bin)
void set_entropy_from_utherm(double utherm, int i)
TimeBinData TimeBinsGravity
void print_particle_info(int i)
Print information relative to a particle / cell to standard output.
void dump_particles(void)
void timebins_get_bin_and_do_validity_checks(integertime ti_step, int *bin_new, int bin_old)
TimeBinData TimeBinsHydro
integertime find_next_sync_point(void)
This function finds the next synchronization point of the system. (i.e. the earliest point of time an...
int TimeBinSynchronized[TIMEBINS]
int getTimeBinSynchronized(int bin)
void make_list_of_active_particles(void)
static bool compare_IDs(const MyIDType &a, const MyIDType &b)
void drift_all_particles(void)
void mark_active_timebins(void)
#define FLAG_BOUNDARYOVERLAP
#define NSOFTCLASSES_HYDRO
#define MAX_DOUBLE_NUMBER
int32_t MySignedIntPosType
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
expr pow(half base, half exp)
void timebins_reallocate(void)
void timebins_allocate(void)
int SofteningClassOfPartType[NTYPES]
int ComovingIntegrationOn
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
double SofteningTable[NSOFTCLASSES+NSOFTCLASSES_HYDRO]
void setSofteningClass(unsigned char softclass)
unsigned char getSofteningClass(void)
unsigned char getTimeBinHydro(void)
unsigned char getType(void)
void myflush(FILE *fstream)