#include <hdf5.h>
#include "allvars.h"
Go to the source code of this file.
Functions | |
void | advance_and_find_timesteps (void) |
void | allocate_commbuffers (void) |
void | allocate_memory (void) |
void | begrun (void) |
int | blockpresent (enum iofields blocknr) |
void | catch_abort (int sig) |
void | catch_fatal (int sig) |
void | check_omega (void) |
void | close_outputfiles (void) |
int | compare_key (const void *a, const void *b) |
void | compute_accelerations (int mode) |
void | compute_global_quantities_of_system (void) |
void | compute_potential (void) |
int | dens_compare_key (const void *a, const void *b) |
void | density (void) |
void | density_decouple (void) |
void | density_evaluate (int i, int mode) |
void | distribute_file (int nfiles, int firstfile, int firsttask, int lasttask, int *filenr, int *master, int *last) |
double | dmax (double, double) |
double | dmin (double, double) |
void | do_box_wrapping (void) |
void | domain_Decomposition (void) |
int | domain_compare_key (const void *a, const void *b) |
int | domain_compare_toplist (const void *a, const void *b) |
void | domain_countToGo (void) |
void | domain_decompose (void) |
void | domain_determineTopTree (void) |
void | domain_exchangeParticles (int partner, int sphflag, int send_count, int recv_count) |
void | domain_findExchangeNumbers (int task, int partner, int sphflag, int *send, int *recv) |
void | domain_findExtent (void) |
int | domain_findSplit (int cpustart, int ncpu, int first, int last) |
void | domain_shiftSplit (void) |
void | domain_sumCost (void) |
void | domain_topsplit (int node, peanokey startkey) |
void | domain_topsplit_local (int node, peanokey startkey) |
double | drift_integ (double a, void *param) |
void | dump_particles (void) |
void | empty_read_buffer (enum iofields blocknr, int offset, int pc, int type) |
void | endrun (int) |
void | energy_statistics (void) |
void | every_timestep_stuff (void) |
void | ewald_corr (double dx, double dy, double dz, double *fper) |
void | ewald_force (int ii, int jj, int kk, double x[3], double force[3]) |
void | ewald_init (void) |
double | ewald_pot_corr (double dx, double dy, double dz) |
double | ewald_psi (double x[3]) |
void | fill_Tab_IO_Labels (void) |
void | fill_write_buffer (enum iofields blocknr, int *pindex, int pc, int type) |
void | find_dt_displacement_constraint (double hfac) |
int | find_files (char *fname) |
int | find_next_outputtime (int time) |
void | find_next_sync_point_and_drift (void) |
void | force_create_empty_nodes (int no, int topnode, int bits, int x, int y, int z, int *nodecount, int *nextfree) |
void | force_exchange_pseudodata (void) |
void | force_flag_localnodes (void) |
void | force_insert_pseudo_particles (void) |
void | force_setupnonrecursive (int no) |
void | force_treeallocate (int maxnodes, int maxpart) |
int | force_treebuild (int npart) |
int | force_treebuild_single (int npart) |
int | force_treeevaluate (int target, int mode, double *ewaldcountsum) |
int | force_treeevaluate_direct (int target, int mode) |
int | force_treeevaluate_ewald_correction (int target, int mode, double pos_x, double pos_y, double pos_z, double aold) |
void | force_treeevaluate_potential (int target, int type) |
void | force_treeevaluate_potential_shortrange (int target, int mode) |
int | force_treeevaluate_shortrange (int target, int mode) |
void | force_treefree (void) |
void | force_treeupdate_pseudos (void) |
void | force_update_hmax (void) |
void | force_update_len (void) |
void | force_update_node (int no, int flag) |
void | force_update_node_hmax_local (void) |
void | force_update_node_hmax_toptree (void) |
void | force_update_node_len_local (void) |
void | force_update_node_len_toptree (void) |
void | force_update_node_recursive (int no, int sib, int father) |
void | force_update_pseudoparticles (void) |
void | force_update_size_of_parent_node (int no) |
void | free_memory (void) |
int | get_bytes_per_blockelement (enum iofields blocknr) |
void | get_dataset_name (enum iofields blocknr, char *buf) |
int | get_datatype_in_block (enum iofields blocknr) |
double | get_drift_factor (int time0, int time1) |
double | get_gravkick_factor (int time0, int time1) |
double | get_hydrokick_factor (int time0, int time1) |
int | get_particles_in_block (enum iofields blocknr, int *typelist) |
double | get_random_number (int id) |
int | get_timestep (int p, double *a, int flag) |
int | get_values_per_blockelement (enum iofields blocknr) |
int | grav_tree_compare_key (const void *a, const void *b) |
void | gravity_forcetest (void) |
void | gravity_tree (void) |
void | gravity_tree_shortrange (void) |
double | gravkick_integ (double a, void *param) |
int | hydro_compare_key (const void *a, const void *b) |
void | hydro_evaluate (int target, int mode) |
void | hydro_force (void) |
double | hydrokick_integ (double a, void *param) |
int | imax (int, int) |
int | imin (int, int) |
void | init (void) |
void | init_drift_table (void) |
void | init_peano_map (void) |
void | long_range_force (void) |
void | long_range_init (void) |
void | long_range_init_regionsize (void) |
void | move_particles (int time0, int time1) |
size_t | my_fread (void *ptr, size_t size, size_t nmemb, FILE *stream) |
size_t | my_fwrite (void *ptr, size_t size, size_t nmemb, FILE *stream) |
int | ngb_clear_buf (float searchcenter[3], float hguess, int numngb) |
void | ngb_treeallocate (int npart) |
void | ngb_treebuild (void) |
int | ngb_treefind_pairs (float searchcenter[3], float hsml, int *startnode) |
int | ngb_treefind_variable (float searchcenter[3], float hguess, int *startnode) |
void | ngb_treefree (void) |
void | ngb_treesearch (int) |
void | ngb_treesearch_pairs (int) |
void | ngb_update_nodes (void) |
void | open_outputfiles (void) |
peanokey | peano_hilbert_key (int x, int y, int z, int bits) |
void | peano_hilbert_order (void) |
void | pm_init_nonperiodic (void) |
void | pm_init_nonperiodic_allocate (int dimprod) |
void | pm_init_nonperiodic_free (void) |
void | pm_init_periodic (void) |
void | pm_init_periodic_allocate (int dimprod) |
void | pm_init_periodic_free (void) |
void | pm_init_regionsize (void) |
void | pm_setup_nonperiodic_kernel (void) |
int | pmforce_nonperiodic (int grnr) |
void | pmforce_periodic (void) |
void | pmpotential_nonperiodic (int grnr) |
void | pmpotential_periodic (void) |
double | pow (double, double) |
void | read_file (char *fname, int readTask, int lastTask) |
void | read_header_attributes_in_hdf5 (char *fname) |
void | read_ic (char *fname) |
int | read_outputlist (char *fname) |
void | read_parameter_file (char *fname) |
void | readjust_timebase (double TimeMax_old, double TimeMax_new) |
void | reorder_gas (void) |
void | reorder_particles (void) |
void | restart (int mod) |
void | run (void) |
void | savepositions (int num) |
double | second (void) |
void | seed_glass (void) |
void | set_random_numbers (void) |
void | set_softenings (void) |
void | set_units (void) |
void | setup_smoothinglengths (void) |
void | statistics (void) |
void | terminate_processes (void) |
double | timediff (double t0, double t1) |
void | write_header_attributes_in_hdf5 (hid_t handle) |
void | write_file (char *fname, int readTask, int lastTask) |
void | write_pid_file (void) |
Definition in file proto.h.
|
This function advances the system in momentum space, i.e. it does apply the 'kick' operation after the forces have been computed. Additionally, it assigns new timesteps to particles. At start-up, a half-timestep is carried out, as well as at the end of the simulation. In between, the half-step kick that ends the previous timestep and the half-step kick for the new timestep are combined into one operation. < adiabatic index of simulated gas < adiabatic index of simulated gas < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < adiabatic index of simulated gas < adiabatic index of simulated gas < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 Definition at line 24 of file timestep.c. References All, global_data_all_processes::ComovingIntegrationOn, global_data_all_processes::CPU_TimeLine, NODE::d, sph_particle_data::DtEntropy, sph_particle_data::Entropy, Extnodes, Father, NODE::father, find_dt_displacement_constraint(), Flag_FullStep, FLOAT, global_data_all_processes::G, GAMMA, GAMMA_MINUS1, get_gravkick_factor(), get_hydrokick_factor(), get_random_number(), get_timestep(), particle_data::GravAccel, particle_data::GravPM, global_data_all_processes::Hubble, sph_particle_data::HydroAccel, particle_data::Mass, NODE::mass, global_data_all_processes::MinEgySpec, Nodes, global_data_all_processes::NumForcesSinceLastDomainDecomp, global_data_all_processes::Omega0, global_data_all_processes::OmegaLambda, P, global_data_all_processes::PM_Ti_begstep, global_data_all_processes::PM_Ti_endstep, particle_data::Pos, pow(), second(), SphP, ThisTask, particle_data::Ti_begstep, global_data_all_processes::Ti_Current, particle_data::Ti_endstep, global_data_all_processes::Time, TIMEBASE, global_data_all_processes::Timebase_interval, timediff(), global_data_all_processes::TotNumPart, global_data_all_processes::TreeDomainUpdateFrequency, particle_data::Type, NODE::u, particle_data::Vel, sph_particle_data::VelPred, and extNODE::vs. Referenced by run(). |
|
Allocates a number of small buffers and arrays, the largest one being the communication buffer. The communication buffer itself is mapped onto various tables used in the different parts of the force algorithms. We further allocate space for the top-level tree nodes, and auxiliary arrays for the domain decomposition algorithm. < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition Definition at line 19 of file allocate.c. References All, global_data_all_processes::BufferSize, global_data_all_processes::BunchSizeDensity, global_data_all_processes::BunchSizeDomain, global_data_all_processes::BunchSizeForce, global_data_all_processes::BunchSizeHydro, CommBuffer, DensDataGet, DensDataIn, DensDataPartialResult, DensDataResult, DomainCount, DomainCountSph, DomainEndList, DomainHmax, DomainKeyBuf, DomainMoment, DomainNodeIndex, DomainPartBuf, DomainSphBuf, DomainStartList, DomainTask, DomainTreeNodeLen, DomainWork, endrun(), Exportflag, FLOAT, GravDataGet, GravDataIn, GravDataIndexTable, GravDataOut, GravDataResult, HydroDataGet, HydroDataIn, HydroDataPartialResult, HydroDataResult, MAXTOPNODES, NTask, peanokey, ThisTask, and TopNodes. Referenced by begrun(). |
|
This routine allocates memory for particle storage, both the collisionless and the SPH particles. Definition at line 103 of file allocate.c. References All, endrun(), global_data_all_processes::MaxPart, global_data_all_processes::MaxPartSph, P, SphP, and ThisTask. Referenced by read_file(), and restart(). |
|
|
This function tells whether or not a given block in the output file is present, depending on the type of simulation run and the compile-time options. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 532 of file io.c. Referenced by read_file(), and write_file(). |
|
|
|
|
|
This routine computes the mass content of the box and compares it to the specified value of Omega-matter. If discrepant, the run is terminated. Definition at line 160 of file init.c. References All, P, and ThisTask. Referenced by init(). |
|
This function closes the global log-files. Definition at line 257 of file begrun.c. References FdCPU, FdEnergy, FdForceTest, FdInfo, FdTimings, and ThisTask. Referenced by run(). |
|
This function is a comparison kernel for sorting the Peano-Hilbert keys. Definition at line 98 of file peano.c. Referenced by peano_hilbert_order(). |
|
This routine computes the accelerations for all active particles. First, the long-range PM force is computed if the TreePM algorithm is used and a "big" PM step is done. Next, the gravitational tree forces are computed. This also constructs the tree, if needed. If gas particles are present, the density-loop for active SPH particles is carried out. This includes an iteration on the correct number of neighbours. Finally, the hydrodynamical forces are added. Definition at line 24 of file accel.c. References All, global_data_all_processes::CPU_Gravity, global_data_all_processes::CPU_Hydro, global_data_all_processes::CPU_PM, global_data_all_processes::CPU_Predict, density(), force_update_hmax(), gravity_forcetest(), gravity_tree(), hydro_force(), long_range_force(), global_data_all_processes::PM_Ti_endstep, second(), ThisTask, global_data_all_processes::Ti_Current, timediff(), global_data_all_processes::TotN_gas, and global_data_all_processes::TypeOfOpeningCriterion. Referenced by run(). |
|
|
|
This routine is a comparison kernel used in a sort routine to group particles that are exported to the same processor. Definition at line 607 of file density.c. Referenced by density(). |
|
|
|
|
|
This function assigns a certain number of files to processors, such that each processor is exactly assigned to one file, and the number of cpus per file is as homogenous as possible. The number of files may at most be equal to the number of processors. Definition at line 725 of file read_ic.c. References ThisTask. Referenced by read_ic(), and savepositions(). |
|
returns the maximum of two double Definition at line 95 of file system.c. Referenced by density(), domain_findSplit(), domain_shiftSplit(), fill_write_buffer(), and read_ic(). |
|
returns the minimum of two double Definition at line 105 of file system.c. Referenced by hydro_evaluate(). |
|
This function makes sure that all particle coordinates (Pos) are periodically mapped onto the interval [0, BoxSize]. After this function has been called, a new domain decomposition should be done, which will also force a new tree construction. Definition at line 103 of file predict.c. References All, global_data_all_processes::BoxSize, P, and particle_data::Pos. Referenced by domain_Decomposition(). |
|
This is a comparison kernel used in a sort routine. |
|
This is a comparison kernel used in a sort routine. |
|
This function determines how many particles that are currently stored on the local CPU have to be moved off according to the domain decomposition. Definition at line 729 of file domain.c. References topnode_data::Daughter, DomainTask, Key, topnode_data::Leaf, NTask, P, topnode_data::Size, topnode_data::StartKey, TopNodes, and particle_data::Type. Referenced by domain_decompose(). |
|
This function carries out the actual domain decomposition for all particle types. It will try to balance the work-load for each domain, as estimated based on the P[i]-GravCost values. The decomposition will respect the maximum allowed memory-imbalance given by the value of PartAllocFactor. Definition at line 154 of file domain.c. References All, domain_countToGo(), domain_determineTopTree(), domain_exchangeParticles(), domain_findExchangeNumbers(), domain_findExtent(), domain_findSplit(), domain_shiftSplit(), domain_sumCost(), DomainEndList, DomainMyLast, DomainMyStart, DomainStartList, endrun(), NTask, NTopleaves, Ntype, NtypeLocal, P, global_data_all_processes::SofteningTable, ThisTask, global_data_all_processes::TotN_gas, and particle_data::Type. Referenced by domain_Decomposition(). |
|
This is the main routine for the domain decomposition. It acts as a driver routine that allocates various temporary buffers, maps the particles back onto the periodic box if needed, and then does the domain decomposition, and a final Peano-Hilbert order of all particles as a tuning measure. Definition at line 62 of file domain.c. References All, global_data_all_processes::CPU_Domain, global_data_all_processes::CPU_Peano, do_box_wrapping(), domain_decompose(), Key, KeySorted, global_data_all_processes::MaxPart, global_data_all_processes::MaxPartSph, N_gas, NTask, global_data_all_processes::NumForcesSinceLastDomainDecomp, NumPart, peano_hilbert_order(), peanokey, global_data_all_processes::PM_Ti_endstep, second(), ThisTask, global_data_all_processes::Ti_Current, timediff(), global_data_all_processes::TotNumPart, global_data_all_processes::TreeDomainUpdateFrequency, and TreeReconstructFlag. |
|
This function constructs the global top-level tree node that is used for the domain decomposition. This is done by considering the string of Peano-Hilbert keys for all particles, which is recursively chopped off in pieces of eight segments until each segment holds at most a certain number of particles. < Bits per dimension available for Peano-Hilbert order. Note: If peanokey is defined as type int, the allowed maximum is 10. If 64-bit integers are used, the maximum is 21 < Bits per dimension available for Peano-Hilbert order. Note: If peanokey is defined as type int, the allowed maximum is 10. If 64-bit integers are used, the maximum is 21 < The number of different Peano-Hilbert cells < Bits per dimension available for Peano-Hilbert order. Note: If peanokey is defined as type int, the allowed maximum is 10. If 64-bit integers are used, the maximum is 21 < The number of different Peano-Hilbert cells Definition at line 898 of file domain.c. References BITS_PER_DIMENSION, DomainCorner, DomainFac, Key, KeySorted, P, and peano_hilbert_key(). Referenced by domain_decompose(). |
|
This function exchanges particles between two CPUs according to the domain split. In doing this, the memory boundaries which may restrict the exhange process are observed. Definition at line 595 of file domain.c. References topnode_data::Daughter, DomainKeyBuf, DomainPartBuf, DomainSphBuf, DomainTask, endrun(), Key, topnode_data::Leaf, N_gas, NumPart, P, peanokey, topnode_data::Size, SphP, topnode_data::StartKey, TAG_KEY, TAG_PDATA, TAG_SPHDATA, ThisTask, TopNodes, and particle_data::Type. Referenced by domain_decompose(). |
|
This function counts how many particles have to be exchanged between two CPUs according to the domain split. If the CPUs are already quite full and hold data from other CPUs as well, not all the particles may be exchanged at once. In this case the communication phase has to be repeated, until enough of the third-party particles have been moved away such that the decomposition can be completed. Definition at line 534 of file domain.c. References All, global_data_all_processes::BunchSizeDomain, imin(), global_data_all_processes::MaxPart, global_data_all_processes::MaxPartSph, and NTask. Referenced by domain_decompose(). |
|
This routine finds the extent of the global domain grid. < Bits per dimension available for Peano-Hilbert order. Note: If peanokey is defined as type int, the allowed maximum is 10. If 64-bit integers are used, the maximum is 21 Definition at line 845 of file domain.c. References DomainCenter, DomainCorner, DomainFac, DomainLen, P, and particle_data::Pos. Referenced by domain_decompose(). |
|
This function tries to find a split point in a range of cells in the domain-grid. The range of cells starts at 'first', and ends at 'last' (inclusively). The number of cpus that holds the range is 'ncpu', with the first cpu given by 'cpustart'. If more than 2 cpus are to be split, the function calls itself recursively. The division tries to achieve a best particle-load balance under the constraint that 'maxload' and 'maxloadsph' may not be exceeded, and that each cpu holds at least one cell from the domaingrid. If such a decomposition cannot be achieved, a non-zero error code is returned. After successful completion, DomainMyStart[] and DomainMyLast[] contain the first and last cell of the domaingrid assigned to the local task for the given type. Also, DomainTask[] contains for each cell the task it was assigned to. Definition at line 327 of file domain.c. References dmax(), DomainCount, DomainCountSph, DomainEndList, DomainStartList, and DomainTask. Referenced by domain_decompose(). |
|
This function tries to improve the domain decomposition found by domain_findSplit() with respect to work-load balance. To this end, the boundaries in the existing domain-split solution (which was found by trying to balance the particle load) are shifted as long as this leads to better work-load while still remaining within the allowed memory-imbalance constraints. Definition at line 447 of file domain.c. References dmax(), DomainCount, DomainCountSph, DomainEndList, DomainStartList, DomainTask, DomainWork, and NTask. Referenced by domain_decompose(). |
|
This routine bins the particles onto the domain-grid, i.e. it sums up the total number of particles and the total amount of work in each of the domain-cells. This information forms the basis for the actual decision on the adopted domain decomposition. Definition at line 786 of file domain.c. References topnode_data::Daughter, domain_walktoptree(), DomainCount, DomainCountSph, DomainWork, particle_data::GravCost, Key, topnode_data::Leaf, NTopleaves, NTopnodes, P, topnode_data::Size, topnode_data::StartKey, ThisTask, particle_data::Ti_begstep, particle_data::Ti_endstep, TopNodes, and particle_data::Type. Referenced by domain_decompose(). |
|
This function is responsible for constructing the global top-level tree segments. Starting from a joint list of all local top-level segments, in which mulitple occurences of the same spatial segment have been combined, a segment is subdivided into 8 pieces recursively until the number of particles in each segment has fallen below All.TotNumPart / (TOPNODEFACTOR * NTask). < Maximum number of nodes in the top-level tree used for domain decomposition Definition at line 1057 of file domain.c. References topnode_data::Blocks, topnode_data::Count, topnode_data::Daughter, endrun(), NTopnodes, topnode_data::Pstart, topnode_data::Size, topnode_data::StartKey, ThisTask, and TopNodes. |
|
This function is responsible for constructing the local top-level Peano-Hilbert segments. A segment is cut into 8 pieces recursively until the number of particles in the segment has fallen below All.TotNumPart / (TOPNODEFACTOR * NTask * NTask). < Maximum number of nodes in the top-level tree used for domain decomposition Definition at line 990 of file domain.c. References topnode_data::Count, topnode_data::Daughter, endrun(), NTopnodes, topnode_data::Pstart, topnode_data::Size, topnode_data::StartKey, ThisTask, and TopNodes. |
|
Integration kernel for drift factor computation. Definition at line 179 of file driftfac.c. References All, global_data_all_processes::Hubble, global_data_all_processes::Omega0, and global_data_all_processes::OmegaLambda. |
|
This function dumps some of the basic particle data to a file. In case the tree construction fails, it is called just before the run terminates with an error message. Examination of the generated file may then give clues to what caused the problem. Definition at line 2763 of file forcetree.c. Referenced by force_create_empty_nodes(), and force_treebuild_single(). |
|
This function reads out the buffer that was filled with particle data, and stores it at the appropriate place in the particle structures. Definition at line 168 of file read_ic.c. References sph_particle_data::Density, sph_particle_data::Entropy, sph_particle_data::Hsml, particle_data::ID, IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, IO_VEL, particle_data::Mass, P, particle_data::Pos, SphP, particle_data::Type, and particle_data::Vel. Referenced by read_file(). |
|
This function aborts the simulations. If a single processors wants an immediate termination, the function needs to be called with ierr>0. A bunch of MPI-error messages may also appear in this case. For ierr=0, MPI is gracefully cleaned up, but this requires that all processors call endrun(). Definition at line 25 of file endrun.c. References terminate_processes(), and ThisTask. Referenced by allocate_commbuffers(), allocate_memory(), density(), domain_decompose(), domain_exchangeParticles(), domain_topsplit(), domain_topsplit_local(), find_files(), force_create_empty_nodes(), force_treeallocate(), force_treebuild_single(), get_particles_in_block(), gravity_forcetest(), init(), long_range_force(), main(), my_fread(), my_fwrite(), ngb_treeallocate(), open_outputfiles(), pm_init_periodic_allocate(), read_file(), read_parameter_file(), readjust_timebase(), restart(), savepositions(), and write_file(). |
|
This routine first calls a computation of various global quantities of the particle distribution, and then writes some statistics about the energies in the various particle components to the file FdEnergy. Definition at line 417 of file run.c. References All, state_of_system::EnergyInt, state_of_system::EnergyIntComp, state_of_system::EnergyKin, state_of_system::EnergyKinComp, state_of_system::EnergyPot, state_of_system::EnergyPotComp, FdEnergy, state_of_system::MassComp, SysState, ThisTask, and global_data_all_processes::Time. Referenced by run(). |
|
|
This function looks up the correction force due to the infinite number of periodic particle/node images. We here use trilinear interpolation to get it from the precomputed tables, which contain one octant around the target particle at the origin. The other octants are obtained from it by exploiting the symmetry properties. Definition at line 2951 of file forcetree.c. |
|
This function computes the force correction term (difference between full force of infinite lattice and nearest image) by Ewald summation. Definition at line 3134 of file forcetree.c. |
|
This function initializes tables with the correction force and the correction potential due to the periodic images of a point mass located at the origin. These corrections are obtained by Ewald summation. (See Hernquist, Bouchet, Suto, ApJS, 1991, 75, 231) The correction fields are used to obtain the full periodic force if periodic boundaries combined with the pure tree algorithm are used. For the TreePM algorithm, the Ewald correction is not used. The correction fields are stored on disk once they are computed. If a corresponding file is found, they are loaded from disk to speed up the initialization. The Ewald summation is done in parallel, i.e. the processors share the work to compute the tables if needed. Definition at line 2802 of file forcetree.c. References ThisTask. Referenced by begrun(). |
|
This function looks up the correction potential due to the infinite number of periodic particle/node images. We here use tri-linear interpolation to get it from the precomputed table, which contains one octant around the target particle at the origin. The other octants are obtained from it by exploiting symmetry properties. Definition at line 3040 of file forcetree.c. References EN. |
|
This function computes the potential correction term by means of Ewald summation. Definition at line 3093 of file forcetree.c. |
|
This function associates a short 4-character block name with each block number. This is stored in front of each block for snapshot FileFormat=2. If one wants to add a new output-block, this function should be augmented accordingly. < total number of defined information blocks for snapshot files. Must be equal to the number of entries in "enum iofields" Definition at line 566 of file io.c. References IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, IO_VEL, iofields, and Tab_IO_Labels. Referenced by read_ic(), and savepositions(). |
|
|
This function computes an upper limit ('dt_displacement') to the global timestep of the system based on the rms velocities of particles. For cosmological simulations, the criterion used is that the rms displacement should be at most a fraction MaxRMSDisplacementFac of the mean particle separation. Note that the latter is estimated using the assigned particle masses, separately for each particle type. If comoving integration is not used, the function imposes no constraint on the timestep.
Definition at line 559 of file timestep.c. Referenced by advance_and_find_timesteps(). |
|
This function determines onto how many files a given snapshot is distributed. Definition at line 616 of file read_ic.c. References All, endrun(), header, global_data_all_processes::ICFormat, io_header::num_files, read_header_attributes_in_hdf5(), and ThisTask. Referenced by read_ic(). |
|
this function returns the next output time that is equal or larger to ti_curr < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 Definition at line 246 of file run.c. References All, global_data_all_processes::ComovingIntegrationOn, global_data_all_processes::OutputListLength, global_data_all_processes::OutputListTimes, global_data_all_processes::Timebase_interval, global_data_all_processes::TimeBegin, and global_data_all_processes::TimeMax. Referenced by begrun(). |
|
This function finds the next synchronization point of the system (i.e. the earliest point of time any of the particles needs a force computation), and drifts the system to this point of time. If the system drifts over the desired time of a snapshot file, the function will drift to this moment, generate an output, and then resume the drift. Definition at line 153 of file run.c. References All, Flag_FullStep, and P. Referenced by run(). |
|
This function recursively creates a set of empty tree nodes which corresponds to the top-level tree for the domain grid. This is done to ensure that this top-level tree is always "complete" so that we can easily associate the pseudo-particles of other CPUs with tree-nodes at a given level in the tree, even when the particle population is so sparse that some of these nodes are actually empty. Definition at line 294 of file forcetree.c. References NODE::center, topnode_data::Daughter, DomainNodeIndex, dump_particles(), endrun(), topnode_data::Leaf, NODE::len, MaxNodes, Nodes, peano_hilbert_key(), NODE::suns, ThisTask, TopNodes, and NODE::u. Referenced by force_treebuild_single(). |
|
This function communicates the values of the multipole moments of the top-level tree-nodes of the domain grid. This data can then be used to update the pseudo-particles on each CPU accordingly. Definition at line 687 of file forcetree.c. References NODE::bitflags, NODE::d, DomainMoment, DomainNodeIndex, Extnodes, NODE::mass, DomainNODE::mass, Nodes, NODE::s, DomainNODE::s, NODE::u, extNODE::vs, and DomainNODE::vs. |
|
This function flags nodes in the top-level tree that are dependent on local particle data. Definition at line 819 of file forcetree.c. References NODE::bitflags, NODE::d, DomainNodeIndex, NODE::father, Nodes, and NODE::u. Referenced by force_treebuild(). |
|
this function inserts pseudo-particles which will represent the mass distribution of the other CPUs. Initially, the mass of the pseudo-particles is set to zero, and their coordinate is set to the center of the domain-cell they correspond to. These quantities will be updated later on. Definition at line 348 of file forcetree.c. References NODE::center, DomainMoment, DomainNodeIndex, DomainNODE::mass, Nodes, and DomainNODE::s. Referenced by force_treebuild_single(). |
|
|
|
This function allocates the memory used for storage of the tree and of auxiliary arrays needed for tree-walk and link-lists. Usually, maxnodes approximately equal to 0.7*maxpart is sufficient to store the tree for up to maxpart particles. < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition Definition at line 2551 of file forcetree.c. References endrun(), MaxNodes, and Nodes_base. Referenced by init(), pmforce_periodic(), pmpotential_periodic(), and restart(). |
|
This function is a driver routine for constructing the gravitational oct-tree, which is done by calling a small number of other functions. Definition at line 61 of file forcetree.c. References All, force_flag_localnodes(), force_treebuild_single(), force_update_pseudoparticles(), Numnodestree, global_data_all_processes::Time, and TimeOfLastTreeConstruction. Referenced by compute_potential(), gravity_tree(), and ngb_treebuild(). |
|
Constructs the gravitational oct-tree. The index convention for accessing tree nodes is the following: the indices 0...NumPart-1 reference single particles, the indices All.MaxPart.... All.MaxPart+nodes-1 reference tree nodes. `Nodes_base' points to the first tree node, while `nodes' is shifted such that nodes[All.MaxPart] gives the first tree node. Finally, node indices with values 'All.MaxPart + MaxNodes' and larger indicate "pseudo particles", i.e. multipole moments of top-level nodes that lie on different CPUs. If such a node needs to be opened, the corresponding particle must be exported to that CPU. The 'Extnodes' structure parallels that of 'Nodes'. Its information is only needed for the SPH part of the computation. (The data is split onto these two structures as a tuning measure. If it is merged into 'Nodes' a somewhat bigger size of the nodes also for gravity would result, which would reduce cache utilization slightly. < Bits per dimension available for Peano-Hilbert order. Note: If peanokey is defined as type int, the allowed maximum is 10. If 64-bit integers are used, the maximum is 21 Definition at line 93 of file forcetree.c. References All, BITS_PER_DIMENSION, NODE::center, NODE::d, topnode_data::Daughter, DomainCenter, DomainCorner, DomainFac, DomainNodeIndex, dump_particles(), endrun(), force_create_empty_nodes(), force_insert_pseudo_particles(), force_update_node_recursive(), global_data_all_processes::ForceSoftening, get_random_number(), particle_data::GravCost, topnode_data::Leaf, NODE::len, MaxNodes, global_data_all_processes::MaxPart, NODE::nextnode, Nextnode, Nodes, P, peano_hilbert_key(), peanokey, particle_data::Pos, topnode_data::Size, topnode_data::StartKey, NODE::suns, ThisTask, TopNodes, particle_data::Type, and NODE::u. Referenced by force_treebuild(). |
|
This routine computes the gravitational force for a given local particle, or for a particle in the communication buffer. Depending on the value of TypeOfOpeningCriterion, either the geometrical BH cell-opening criterion, or the `relative' opening criterion is used. Definition at line 1111 of file forcetree.c. References All, global_data_all_processes::ErrTolForceAcc, sph_particle_data::Hsml, particle_data::OldAcc, P, particle_data::Pos, SphP, and particle_data::Type. Referenced by gravity_tree(). |
|
This function does the force computation with direct summation for the specified particle in the communication buffer. This can be useful for debugging purposes, in particular for explicit checks of the force accuracy. Definition at line 2632 of file forcetree.c. References All, P, particle_data::Pos, and particle_data::Type. Referenced by gravity_forcetest(). |
|
This function computes the Ewald correction, and is needed if periodic boundary conditions together with a pure tree algorithm are used. Note that the ordinary tree walk does not carry out this correction directly as it was done in Gadget-1.1. Instead, the tree is walked a second time. This is actually faster because the "Ewald-Treewalk" can use a different opening criterion than the normal tree walk. In particular, the Ewald correction is negligible for particles that are very close, but it is large for particles that are far away (this is quite different for the normal direct force). So we can here use a different opening criterion. Sufficient accuracy is usually obtained if the node length has dropped to a certain fraction ~< 0.25 of the BoxLength. However, we may only short-cut the interaction list of the normal full Ewald tree walk if we are sure that the whole node and all daughter nodes "lie on the same side" of the periodic boundary, i.e. that the real tree walk would not find a daughter node or particle that was mapped to a different nearest neighbour position when the tree walk would be further refined. Definition at line 1734 of file forcetree.c. References gravdata_in::Acc, All, NODE::bitflags, global_data_all_processes::BoxSize, NODE::center, NODE::d, DomainTask, EN, global_data_all_processes::ErrTolTheta, Exportflag, particle_data::GravAccel, particle_data::GravCost, GravDataResult, NODE::len, NODE::mass, particle_data::Mass, global_data_all_processes::MaxPart, NEAREST, NODE::nextnode, Nextnode, gravdata_in::Ninteractions, Nodes, P, particle_data::Pos, NODE::s, NODE::sibling, gravdata_in::u, NODE::u, and gravdata_in::w. |
|
This routine computes the gravitational potential by walking the tree. The same opening criteria is used as for the gravitational force walk. Definition at line 2000 of file forcetree.c. References All, global_data_all_processes::ErrTolForceAcc, sph_particle_data::Hsml, particle_data::OldAcc, P, particle_data::Pos, SphP, and particle_data::Type. Referenced by compute_potential(). |
|
This function computes the short-range potential when the TreePM algorithm is used. This potential is the Newtonian potential, modified by a complementary error function. Definition at line 2247 of file forcetree.c. References All, global_data_all_processes::ErrTolForceAcc, sph_particle_data::Hsml, particle_data::OldAcc, P, particle_data::Pos, SphP, and particle_data::Type. Referenced by compute_potential(). |
|
In the TreePM algorithm, the tree is walked only locally around the target coordinate. Tree nodes that fall outside a box of half side-length Rcut= RCUT*ASMTH*MeshSize can be discarded. The short-range potential is modified by a complementary error function, multiplied with the Newtonian form. The resulting short-range suppression compared to the Newtonian force is tabulated, because looking up from this table is faster than recomputing the corresponding factor, despite the memory-access panelty (which reduces cache performance) incurred by the table. Definition at line 1393 of file forcetree.c. References All, global_data_all_processes::ErrTolForceAcc, sph_particle_data::Hsml, particle_data::OldAcc, P, particle_data::Pos, SphP, and particle_data::Type. Referenced by gravity_tree(). |
|
This function frees the memory allocated for the tree, i.e. it frees the space allocated by the function force_treeallocate(). Definition at line 2615 of file forcetree.c. Referenced by pmforce_periodic(), and pmpotential_periodic(). |
|
This function updates the top-level tree after the multipole moments of the pseudo-particles have been updated. Definition at line 731 of file forcetree.c. References All, NODE::bitflags, NODE::d, DomainMoment, DomainNodeIndex, Extnodes, NODE::father, FLOAT, global_data_all_processes::ForceSoftening, DomainNODE::mass, NODE::mass, Nodes, DomainNODE::s, NODE::s, NODE::u, DomainNODE::vs, and extNODE::vs. |
|
This function updates the hmax-values in tree nodes that hold SPH particles. These values are needed to find all neighbors in the hydro-force computation. Since the Hsml-values are potentially changed in the SPH-denity computation, force_update_hmax() should be carried out just before the hydrodynamical SPH forces are computed, i.e. after density(). Definition at line 999 of file forcetree.c. References DomainHmax, DomainNodeIndex, and Extnodes. Referenced by compute_accelerations(). |
|
This function updates the side-length of tree nodes in case the tree is not reconstructed, but only drifted. The grouping of particles to tree nodes is not changed in this case, but some tree nodes may need to be enlarged because particles moved out of their original bounds. Definition at line 870 of file forcetree.c. References DomainNodeIndex, DomainTreeNodeLen, and Nodes. Referenced by move_particles(). |
|
|
|
This routine updates the hmax-values of local tree nodes. Definition at line 1036 of file forcetree.c. References NODE::d, Extnodes, NODE::father, Father, extNODE::hmax, sph_particle_data::Hsml, Nodes, SphP, and NODE::u. |
|
This function recursively sets the hmax-values of the top-level tree. Definition at line 1072 of file forcetree.c. References NODE::d, DomainHmax, DomainNodeIndex, Extnodes, NODE::father, extNODE::hmax, Nodes, and NODE::u. |
|
This function recursively enlarges nodes such that they always contain all their daughter nodes and daughter particles. Definition at line 908 of file forcetree.c. References NODE::center, NODE::d, NODE::father, Father, FLOAT, NODE::len, Nodes, P, particle_data::Pos, and NODE::u. |
|
This function recursively enlarges nodes of the top-level tree such that they always contain all their daughter nodes. Definition at line 955 of file forcetree.c. References NODE::center, NODE::d, DomainNodeIndex, DomainTreeNodeLen, NODE::father, FLOAT, NODE::len, Nodes, and NODE::u. |
|
this routine determines the multipole moments for a given internal node and all its subnodes using a recursive computation. The result is stored in the Nodes[] structure in the sequence of this tree-walk. Note that the bitflags-variable for each node is used to store in the lowest bits some special information: Bit 0 flags whether the node belongs to the top-level tree corresponding to the domain decomposition, while Bit 1 signals whether the top-level node is dependent on local mass. If UNEQUALSOFTENINGS is set, bits 2-4 give the particle type with the maximum softening among the particles in the node, and bit 5 flags whether the node contains any particles with lower softening than that. Definition at line 425 of file forcetree.c. References All, NODE::bitflags, NODE::center, NODE::d, Extnodes, NODE::father, FLOAT, global_data_all_processes::ForceSoftening, extNODE::hmax, sph_particle_data::Hsml, particle_data::Mass, NODE::mass, global_data_all_processes::MaxPart, NODE::nextnode, Nextnode, Nodes, P, particle_data::Pos, NODE::s, NODE::sibling, SphP, NODE::suns, particle_data::Type, NODE::u, particle_data::Vel, and extNODE::vs. Referenced by force_treebuild_single(). |
|
This function updates the multipole moments of the pseudo-particles that represent the mass distribution on different CPUs. For that purpose, it first exchanges the necessary data, and then updates the top-level tree accordingly. The detailed implementation of these two tasks is done in separate functions. Definition at line 674 of file forcetree.c. Referenced by force_treebuild(), and move_particles(). |
|
|
|
This routine frees the memory for the particle storage. Note: We don't actually bother to call it in the code... When the program terminats, the memory will be automatically freed by the operating system. Definition at line 144 of file allocate.c. References All, global_data_all_processes::MaxPart, global_data_all_processes::MaxPartSph, P, and SphP. |
|
This function tells the size of one data entry in each of the blocks defined for the output file. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 366 of file io.c. References IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, and IO_VEL. Referenced by read_file(), and write_file(). |
|
This function returns a descriptive character string that describes the name of the block when the HDF5 file format is used. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 614 of file io.c. References IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, and IO_VEL. Referenced by read_file(), and write_file(). |
|
This function returns the type of the data contained in a given block of the output file. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 405 of file io.c. References IO_ID. Referenced by read_file(), and write_file(). |
|
This function integrates the cosmological prefactor for a drift step between time0 and time1. The value returned is *
< length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors Definition at line 67 of file driftfac.c. References All, DRIFT_TABLE_LENGTH, DriftTable, and global_data_all_processes::Timebase_interval. Referenced by move_particles(). |
|
This function integrates the cosmological prefactor for a kick step of the gravitational force. < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors Definition at line 105 of file driftfac.c. References All, DRIFT_TABLE_LENGTH, GravKickTable, and global_data_all_processes::Timebase_interval. Referenced by advance_and_find_timesteps(), compute_global_quantities_of_system(), fill_write_buffer(), and move_particles(). |
|
This function integrates the cosmological prefactor for a kick step of the hydrodynamical force. < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors Definition at line 142 of file driftfac.c. References All, DRIFT_TABLE_LENGTH, HydroKickTable, and global_data_all_processes::Timebase_interval. Referenced by advance_and_find_timesteps(), compute_global_quantities_of_system(), fill_write_buffer(), and move_particles(). |
|
This function determines how many particles there are in a given block, based on the information in the header-structure. It also flags particle types that are present in the block in the typelist array. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 465 of file io.c. References All, endrun(), header, IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, IO_VEL, global_data_all_processes::MassTable, and io_header::npart. Referenced by read_file(), and write_file(). |
|
This routine returns a random number taken from a table of random numbers, which is refilled every timestep. This method is used to allow random number application to particles independent of the number of processors used, and independent of the particular order the particles have. In order to work properly, the particle IDs should be set properly to unique integer values. < gives the length of a table with random numbers, refreshed at every timestep. This is used to allow application of random numbers to a specific particle in a way that is independent of the number of processors used. Definition at line 29 of file system.c. References RndTable. Referenced by advance_and_find_timesteps(), force_treebuild_single(), and gravity_forcetest(). |
|
This function normally (for flag==0) returns the maximum allowed timestep of a particle, expressed in terms of the integer mapping that is used to represent the total simulated timespan. The physical acceleration is returned in `aphys'. The latter is used in conjunction with the PSEUDOSYMMETRIC integration option, which also makes of the second function of get_timestep. When it is called with a finite timestep for flag, it returns the physical acceleration that would lead to this timestep, assuming timestep criterion 0. < adiabatic index of simulated gas < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29
Definition at line 414 of file timestep.c. References particle_data::GravAccel, particle_data::GravPM, sph_particle_data::HydroAccel, P, SphP, and particle_data::Type. Referenced by advance_and_find_timesteps(). |
|
This function informs about the number of elements stored per particle for the given block of the output file. If one wants to add a new output-block, this function should be augmented accordingly. Definition at line 432 of file io.c. References IO_ACCEL, IO_DTENTR, IO_HSML, IO_ID, IO_MASS, IO_POS, IO_POT, IO_RHO, IO_TSTP, IO_U, and IO_VEL. Referenced by read_file(), and write_file(). |
|
This function is used as a comparison kernel in a sort routine. It is used to group particles in the communication buffer that are going to be sent to the same CPU. Definition at line 511 of file gravtree.c. Referenced by compute_potential(), gravity_forcetest(), and gravity_tree(). |
|
|
|
|
|
Integration kernel for gravitational kick factor computation. Definition at line 191 of file driftfac.c. References All, global_data_all_processes::Hubble, global_data_all_processes::Omega0, and global_data_all_processes::OmegaLambda. |
|
This is a comparison kernel for a sort routine, which is used to group particles that are going to be exported to the same CPU. Definition at line 563 of file hydra.c. Referenced by hydro_force(). |
|
|
|
Integration kernel for hydrodynamical kick factor computation. < adiabatic index of simulated gas Definition at line 204 of file driftfac.c. References All, GAMMA_MINUS1, global_data_all_processes::Hubble, global_data_all_processes::Omega0, global_data_all_processes::OmegaLambda, and pow(). |
|
returns the maximum of two integers Definition at line 115 of file system.c. Referenced by hydro_evaluate(), and pm_init_periodic_allocate(). |
|
returns the minimum of two integers Definition at line 125 of file system.c. Referenced by domain_findExchangeNumbers(). |
|
|
This function computes look-up tables for factors needed in cosmological integrations. The (simple) integrations are carried out with the GSL library. Separate factors are computed for the "drift", and the gravitational and hydrodynamical "kicks". The lookup-table is used for reasons of speed. < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors < length of the lookup table used to hold the drift and kick factors Definition at line 26 of file driftfac.c. References All, DRIFT_TABLE_LENGTH, DriftTable, GravKickTable, global_data_all_processes::Hubble, HydroKickTable, global_data_all_processes::TimeBegin, global_data_all_processes::TimeMax, and WORKSIZE. Referenced by begrun(). |
|
|
|
This function is a driver routine for the long-range PM force computation. It calls periodic and/or non-periodic FFT routines as needed for the present simulation set-up. Definition at line 56 of file longrange.c. References All, global_data_all_processes::ComovingIntegrationOn, endrun(), particle_data::GravPM, global_data_all_processes::Hubble, global_data_all_processes::Omega0, global_data_all_processes::OmegaLambda, P, pm_init_regionsize(), pm_setup_nonperiodic_kernel(), pmforce_nonperiodic(), pmforce_periodic(), and particle_data::Pos. Referenced by compute_accelerations(). |
|
Calls initializiation routines of periodic or/and non-periodic FFT routines. Definition at line 20 of file longrange.c. References pm_init_nonperiodic(), and pm_init_periodic(). Referenced by begrun(). |
|
This function calls subroutines that determine the spatial region covered by the PM mesh. Definition at line 36 of file longrange.c. References pm_init_regionsize(), pm_setup_nonperiodic_kernel(), and RestartFlag. Referenced by begrun(). |
|
|
This catches I/O errors occuring for fread(). In this case we better stop. Definition at line 1141 of file io.c. References endrun(), and ThisTask. Referenced by read_file(). |
|
This catches I/O errors occuring for my_fwrite(). In this case we better stop. Definition at line 1124 of file io.c. References endrun(), and ThisTask. Referenced by write_file(). |
|
The buffer for the neighbour list has a finite length MAX_NGB. For a large search region, this buffer can get full, in which case this routine can be called to eliminate some of the superfluous particles in the "corners" of the search box - only the ones in the inscribed sphere need to be kept. Definition at line 320 of file ngb.c. References FLOAT, NGB_PERIODIC_X, NGB_PERIODIC_Y, NGB_PERIODIC_Z, Ngblist, P, and particle_data::Pos. Referenced by ngb_treefind_variable(). |
|
Allocates memory for the neighbour list buffer. Definition at line 358 of file ngb.c. References All, boxHalf_X, boxHalf_Y, boxHalf_Z, global_data_all_processes::BoxSize, boxSize_X, boxSize_Y, boxSize_Z, endrun(), Ngblist, and ThisTask. |
|
This function constructs the neighbour tree. To this end, we actually need to construct the gravitational tree, because we use it now for the neighbour search. Definition at line 403 of file ngb.c. References force_treebuild(), N_gas, and ThisTask. Referenced by init(). |
|
This routine finds all neighbours `j' that can interact with the particle `i' in the communication buffer. Note that an interaction can take place if OR if . In the range-search this is taken into account, i.e. it is guaranteed that all particles are found that fulfil this condition, including the (more difficult) second part of it. For this purpose, each node knows the maximum h occuring among the particles it represents. < defines maximum length of neighbour list Definition at line 64 of file ngb.c. References All, DomainTask, Exportflag, Extnodes, FLOAT, extNODE::hmax, sph_particle_data::Hsml, global_data_all_processes::MaxPart, Nextnode, NGB_PERIODIC_X, NGB_PERIODIC_Y, NGB_PERIODIC_Z, Ngblist, Nodes, P, particle_data::Pos, SphP, and particle_data::Type. Referenced by hydro_evaluate(). |
|
This function returns neighbours with distance <= hsml and returns them in Ngblist. Actually, particles in a box of half side length hsml are returned, i.e. the reduction to a sphere still needs to be done in the calling routine. < defines maximum length of neighbour list < defines maximum length of neighbour list Definition at line 194 of file ngb.c. References All, NODE::d, DomainTask, Exportflag, FLOAT, global_data_all_processes::MaxPart, NODE::nextnode, Nextnode, ngb_clear_buf(), NGB_PERIODIC_X, NGB_PERIODIC_Y, NGB_PERIODIC_Z, Ngblist, Nodes, P, particle_data::Pos, NODE::sibling, ThisTask, and particle_data::Type. Referenced by density_evaluate(). |
|
free memory allocated for neighbour list buffer. Definition at line 394 of file ngb.c. References Ngblist. |
|
|
|
|
|
|
|
This function opens various log-files that report on the status and performance of the simulstion. On restart from restart-files (start-option 1), the code will append to these files. Definition at line 199 of file begrun.c. References All, global_data_all_processes::CpuFile, endrun(), global_data_all_processes::EnergyFile, FdCPU, FdEnergy, FdForceTest, FdInfo, FdTimings, global_data_all_processes::InfoFile, global_data_all_processes::OutputDir, RestartFlag, ThisTask, and global_data_all_processes::TimingsFile. Referenced by begrun(). |
|
This function computes a Peano-Hilbert key for an integer triplet (x,y,z), with x,y,z in the range between 0 and 2^bits-1. Definition at line 263 of file peano.c. References peanokey. Referenced by domain_determineTopTree(), force_create_empty_nodes(), and force_treebuild_single(). |
|
This function puts the particles into Peano-Hilbert order by sorting them according to their keys. The latter half already been computed in the domain decomposition. Since gas particles need to stay at the beginning of the particle list, they are sorted as a separate block. Definition at line 34 of file peano.c. References compare_key(), Key, N_gas, NumPart, reorder_gas(), reorder_particles(), and ThisTask. Referenced by domain_Decomposition(). |
|
Referenced by long_range_init(). |
|
|
|
|
|
This routines generates the FFTW-plans to carry out the parallel FFTs later on. Some auxiliary variables are also initialized. < @-ASMTH gives the scale of the short-range/long-range force split in units of FFT-mesh cells < @-RCUT gives the maximum distance (in units of the scale used for the force split) out to which short-range forces are evaluated in the short-range tree walk. Definition at line 55 of file pm_periodic.c. References All, ASMTH, global_data_all_processes::Asmth, global_data_all_processes::BoxSize, NTask, RCUT, global_data_all_processes::Rcut, and ThisTask. Referenced by long_range_init(). |
|
This function allocates the memory neeed to compute the long-range PM force. Three fields are used, one to hold the density (and its FFT, and then the real-space potential), one to hold the force field obtained by finite differencing, and finally a workspace field, which is used both as workspace for the parallel FFT, and as buffer for the communication algorithm used in the force computation. Definition at line 114 of file pm_periodic.c. References endrun(), imax(), and ThisTask. Referenced by pmforce_periodic(), and pmpotential_periodic(). |
|
This routine frees the space allocated for the parallel FFT algorithm. Definition at line 161 of file pm_periodic.c. Referenced by pmforce_periodic(), and pmpotential_periodic(). |
|
Referenced by long_range_force(), and long_range_init_regionsize(). |
|
Referenced by long_range_force(), and long_range_init_regionsize(). |
|
Referenced by long_range_force(). |
|
Calculates the long-range periodic force given the particle positions using the PM method. The force is Gaussian filtered with Asmth, given in mesh-cell units. We carry out a CIC charge assignment, and compute the potenial by Fourier transform methods. The potential is finite differenced using a 4-point finite differencing formula, and the forces are interpolated tri-linearly to the particle positions. The CIC kernel is deconvolved. Note that the particle distribution is not in the slab decomposition that is used for the FFT. Instead, overlapping patches between local domains and FFT slabs are communicated as needed. Definition at line 181 of file pm_periodic.c. References All, global_data_all_processes::Asmth, global_data_all_processes::BoxSize, force_treeallocate(), force_treefree(), global_data_all_processes::G, particle_data::GravPM, particle_data::Mass, global_data_all_processes::MaxPart, global_data_all_processes::NumForcesSinceLastDomainDecomp, P, pm_init_periodic_allocate(), pm_init_periodic_free(), PMGRID2, particle_data::Pos, TAG_PERIODIC_A, TAG_PERIODIC_B, ThisTask, global_data_all_processes::TotNumPart, global_data_all_processes::TreeAllocFactor, and global_data_all_processes::TreeDomainUpdateFrequency. Referenced by long_range_force(). |
|
Referenced by compute_potential(). |
|
Calculates the long-range potential using the PM method. The potential is Gaussian filtered with Asmth, given in mesh-cell units. We carry out a CIC charge assignment, and compute the potenial by Fourier transform methods. The CIC kernel is deconvolved. Definition at line 694 of file pm_periodic.c. References All, global_data_all_processes::Asmth, global_data_all_processes::BoxSize, force_treeallocate(), force_treefree(), global_data_all_processes::G, particle_data::Mass, global_data_all_processes::MaxPart, global_data_all_processes::NumForcesSinceLastDomainDecomp, P, pm_init_periodic_allocate(), pm_init_periodic_free(), PMGRID2, particle_data::Pos, particle_data::Potential, TAG_PERIODIC_C, TAG_PERIODIC_D, ThisTask, global_data_all_processes::TotNumPart, global_data_all_processes::TreeAllocFactor, and global_data_all_processes::TreeDomainUpdateFrequency. Referenced by compute_potential(). |
|
|
|
This function reads the header information in case the HDF5 file format is used. Definition at line 765 of file read_ic.c. References io_header::flag_entropy_instead_u, header, io_header::mass, io_header::npart, io_header::npartTotal, io_header::npartTotalHighWord, and io_header::num_files. Referenced by find_files(), and read_file(). |
|
This function reads initial conditions, in one of the three possible file formats currently supported by Gadget. Note: When a snapshot file is started from initial conditions (start-option 0), not all the information in the header is used, in particular, the STARTING TIME needs to be set in the parameterfile. Also, for gas particles, only the internal energy is read, the density and mean molecular weight will be recomputed by the code. When InitGasTemp>0 is given, the gas temperature will be initialzed to this value assuming a mean colecular weight either corresponding to complete neutrality, or full ionization. However, when the code is started with start-option 2, then all the this data in the snapshot files is preserved, i.e. this is also the way to resume a simulation from a snapshot file in case a regular restart file is not available. < adiabatic index of simulated gas < mass fraction of hydrogen, relevant only for radiative cooling < mass fraction of hydrogen, relevant only for radiative cooling Definition at line 31 of file read_ic.c. References All, BOLTZMANN, distribute_file(), dmax(), sph_particle_data::Entropy, fill_Tab_IO_Labels(), find_files(), global_data_all_processes::ICFormat, global_data_all_processes::InitGasTemp, particle_data::Mass, global_data_all_processes::MassTable, global_data_all_processes::MinEgySpec, N_gas, NTask, global_data_all_processes::NumFilesWrittenInParallel, NumPart, P, read_file(), RestartFlag, SphP, ThisTask, global_data_all_processes::TotNumPart, particle_data::Type, global_data_all_processes::UnitEnergy_in_cgs, and global_data_all_processes::UnitMass_in_g. Referenced by init(). |
|
this function reads a table with a list of desired output times. The table does not have to be ordered in any way, but may not contain more than MAXLEN_OUTPUTLIST entries. < maxmimum number of entries in list of snapshot output times Definition at line 759 of file begrun.c. References All, global_data_all_processes::OutputListLength, and global_data_all_processes::OutputListTimes. Referenced by read_parameter_file(). |
|
|
If a restart from restart-files is carried out where the TimeMax variable is increased, then the integer timeline needs to be adjusted. The approach taken here is to reduce the resolution of the integer timeline by factors of 2 until the new final time can be reached within TIMEBASE. < The simulated timespan is mapped onto the integer interval [0,TIMESPAN], where TIMESPAN needs to be a power of 2. Note that (1<<28) corresponds to 2^29 Definition at line 793 of file begrun.c. References All, global_data_all_processes::ComovingIntegrationOn, endrun(), P, global_data_all_processes::PM_Ti_begstep, global_data_all_processes::PM_Ti_endstep, ThisTask, particle_data::Ti_begstep, global_data_all_processes::Ti_Current, particle_data::Ti_endstep, global_data_all_processes::Timebase_interval, global_data_all_processes::TimeBegin, and global_data_all_processes::TimeMax. Referenced by begrun(). |
|
This function brings the gas particles into the same order as the sorted keys. (The sort is first done only on the keys themselves and done directly on the gas particles in order to reduce the amount of data that needs to be moved in memory. Only once the order is established, the gas particles are rearranged, such that each particle has to be moved at most once.) Definition at line 117 of file peano.c. Referenced by peano_hilbert_order(). |
|
This function brings the collisionless particles into the same order as the sorted keys. (The sort is first done only on the keys themselves and done directly on the particles in order to reduce the amount of data that needs to be moved in memory. Only once the order is established, the particles are rearranged, such that each particle has to be moved at most once.) Definition at line 166 of file peano.c. References P. Referenced by peano_hilbert_order(). |
|
This function reads or writes the restart files. Each processor writes its own restart file, with the I/O being done in parallel. To avoid congestion of the disks you can tell the program to restrict the number of files that are simultaneously written to NumFilesWrittenInParallel. If modus>0 the restart()-routine reads, if modus==0 it writes a restart file. < defines maximum length of neighbour list < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition < Maximum number of nodes in the top-level tree used for domain decomposition Definition at line 35 of file restart.c. References All, allocate_memory(), NODE::d, DomainCenter, DomainCorner, DomainEndList, DomainFac, DomainHmax, DomainLen, DomainMoment, DomainMyLast, DomainMyStart, DomainNodeIndex, DomainStartList, DomainTask, DomainTreeNodeLen, endrun(), Extnodes_base, Father, NODE::father, FLOAT, force_treeallocate(), MAX_NGB, MaxNodes, global_data_all_processes::MaxPart, global_data_all_processes::MaxPartSph, MAXTOPNODES, N_gas, Nextnode, NODE::nextnode, ngb_treeallocate(), Nodes_base, NTask, global_data_all_processes::NumFilesWrittenInParallel, Numnodestree, NumPart, global_data_all_processes::OutputDir, P, global_data_all_processes::PartAllocFactor, random_generator, global_data_all_processes::RestartFile, NODE::sibling, SphP, ThisTask, global_data_all_processes::Time, global_data_all_processes::TotN_gas, global_data_all_processes::TotNumPart, global_data_all_processes::TreeAllocFactor, and NODE::u. |
|
|
This function writes a snapshot of the particle distribution to one or several files using the selected file format. If NumFilesPerSnapshot>1, the snapshot is distributed onto several files, several of them can be written simultaneously (up to NumFilesWrittenInParallel). Each file contains data from a group of processors. Definition at line 33 of file io.c. References All, global_data_all_processes::CPU_Snapshot, distribute_file(), endrun(), fill_Tab_IO_Labels(), NTask, global_data_all_processes::NumFilesPerSnapshot, global_data_all_processes::NumFilesWrittenInParallel, global_data_all_processes::NumForcesSinceLastDomainDecomp, global_data_all_processes::OutputDir, P, second(), global_data_all_processes::SnapFormat, global_data_all_processes::SnapshotFileBase, ThisTask, timediff(), global_data_all_processes::TotNumPart, global_data_all_processes::TreeDomainUpdateFrequency, particle_data::Type, and write_file(). Referenced by run(). |
|
returns the number of cpu-ticks in seconds that have elapsed, or the wall-clock time obtained with MPI_Wtime(). Definition at line 53 of file system.c. Referenced by advance_and_find_timesteps(), compute_accelerations(), compute_potential(), density(), domain_Decomposition(), gravity_forcetest(), gravity_tree(), hydro_force(), main(), move_particles(), run(), and savepositions(). |
|
Referenced by init(). |
|
This routine fills the random number table. < gives the length of a table with random numbers, refreshed at every timestep. This is used to allow application of random numbers to a specific particle in a way that is independent of the number of processors used. Definition at line 39 of file system.c. References RndTable. Referenced by begrun(). |
|
|
Computes conversion factors between internal code units and the cgs-system. < Gravitational constant (in cgs units) < mass fraction of hydrogen, relevant only for radiative cooling < adiabatic index of simulated gas Definition at line 149 of file begrun.c. References All, BOLTZMANN, global_data_all_processes::G, GRAVITY, global_data_all_processes::GravityConstantInternal, HUBBLE, global_data_all_processes::Hubble, global_data_all_processes::MinEgySpec, global_data_all_processes::MinGasTemp, pow(), ThisTask, global_data_all_processes::UnitCoolingRate_in_cgs, global_data_all_processes::UnitDensity_in_cgs, global_data_all_processes::UnitEnergy_in_cgs, global_data_all_processes::UnitLength_in_cm, global_data_all_processes::UnitMass_in_g, global_data_all_processes::UnitPressure_in_cgs, global_data_all_processes::UnitTime_in_Megayears, global_data_all_processes::UnitTime_in_s, and global_data_all_processes::UnitVelocity_in_cm_per_s. Referenced by begrun(). |
|
This function is used to find an initial smoothing length for each SPH particle. It guarantees that the number of neighbours will be between desired_ngb-MAXDEV and desired_ngb+MAXDEV. For simplicity, a first guess of the smoothing length is provided to the function density(), which will then iterate if needed to find the right smoothing length. Definition at line 197 of file init.c. References All, NODE::d, global_data_all_processes::DesNumNgb, NODE::father, Father, sph_particle_data::Hsml, NODE::len, NODE::mass, particle_data::Mass, Nodes, P, pow(), RestartFlag, SphP, and NODE::u. Referenced by init(). |
|
|
|
Referenced by endrun(). |
|
returns the time difference between two measurements obtained with second(). The routine takes care of the possible overflow of the tick counter on 32bit systems, but depending on the system, this may not always work properly. Similarly, in some MPI implementations, the MPI_Wtime() function may also overflow, in which case a negative time difference would be returned. The routine returns instead a time difference equal to 0. Definition at line 74 of file system.c. References pow(). Referenced by advance_and_find_timesteps(), compute_accelerations(), compute_potential(), density(), domain_Decomposition(), gravity_forcetest(), gravity_tree(), hydro_force(), main(), move_particles(), run(), and savepositions(). |
|
This function writes an actual snapshot file containing the data from processors 'writeTask' to 'lastTask'. 'writeTask' is the one that actually writes. Each snapshot file contains a header first, then particle positions, velocities and ID's. Particle masses are written only for those particle types with zero entry in MassTable. After that, first the internal energies u, and then the density is written for the SPH particles. If cooling is enabled, mean molecular weight and neutral hydrogen abundance are written for the gas particles. This is followed by the SPH smoothing length and further blocks of information, depending on included physics and compile-time flags. If HDF5 is used, the header is stored in a group called "/Header", and the particle data is stored separately for each particle type in groups calles "/PartType0", "/PartType1", etc. The sequence of the blocks is unimportant in this case. < Various tags used for labelling MPI messages < Various tags used for labelling MPI messages < total number of defined information blocks for snapshot files. Must be equal to the number of entries in "enum iofields" Definition at line 673 of file io.c. References All, blockpresent(), global_data_all_processes::BoxSize, io_header::BoxSize, global_data_all_processes::BufferSize, CommBuffer, global_data_all_processes::ComovingIntegrationOn, endrun(), fill_write_buffer(), io_header::flag_cooling, io_header::flag_feedback, io_header::flag_metals, io_header::flag_sfr, io_header::flag_stellarage, get_bytes_per_blockelement(), get_dataset_name(), get_datatype_in_block(), get_particles_in_block(), get_values_per_blockelement(), header, global_data_all_processes::HubbleParam, io_header::HubbleParam, iofields, io_header::mass, global_data_all_processes::MassTable, my_fwrite(), io_header::npart, io_header::npartTotal, io_header::npartTotalHighWord, io_header::num_files, global_data_all_processes::NumFilesPerSnapshot, global_data_all_processes::Omega0, io_header::Omega0, global_data_all_processes::OmegaLambda, io_header::OmegaLambda, io_header::redshift, global_data_all_processes::SnapFormat, Tab_IO_Labels, TAG_LOCALN, TAG_N, TAG_NFORTHISTASK, TAG_PDATA, ThisTask, global_data_all_processes::Time, io_header::time, and write_header_attributes_in_hdf5(). Referenced by savepositions(). |
|
This function writes the header information in case HDF5 is selected as file format. Definition at line 1000 of file io.c. Referenced by write_file(). |
|
|