Main Page | Alphabetical List | Data Structures | File List | Data Fields | Globals | Related Pages

global_data_all_processes Struct Reference

#include <allvars.h>


Data Fields

long long TotNumPart
long long TotN_gas
int MaxPart
int MaxPartSph
double BoxSize
int ICFormat
int SnapFormat
int NumFilesPerSnapshot
int NumFilesWrittenInParallel
int BufferSize
int BunchSizeForce
int BunchSizeDensity
int BunchSizeHydro
int BunchSizeDomain
double PartAllocFactor
double TreeAllocFactor
double DesNumNgb
double MaxNumNgbDeviation
double ArtBulkViscConst
double InitGasTemp
double MinGasTemp
double MinEgySpec
long long TotNumOfForces
long long NumForcesSinceLastDomainDecomp
double G
double UnitTime_in_s
double UnitMass_in_g
double UnitVelocity_in_cm_per_s
double UnitLength_in_cm
double UnitPressure_in_cgs
double UnitDensity_in_cgs
double UnitCoolingRate_in_cgs
double UnitEnergy_in_cgs
double UnitTime_in_Megayears
double GravityConstantInternal
double Hubble
double Omega0
double OmegaLambda
double OmegaBaryon
double HubbleParam
int ComovingIntegrationOn
int PeriodicBoundariesOn
int ResubmitOn
int TypeOfOpeningCriterion
int TypeOfTimestepCriterion
int OutputListOn
int SnapshotFileCount
double TimeBetSnapshot
double TimeOfFirstSnapshot
double CpuTimeBetRestartFile
double TimeLastRestartFile
double TimeBetStatistics
double TimeLastStatistics
int NumCurrentTiStep
double Time
double TimeBegin
double TimeStep
double TimeMax
double Timebase_interval
int Ti_Current
int Ti_nextoutput
int PM_Ti_endstep
int PM_Ti_begstep
double Asmth [2]
double Rcut [2]
double Corner [2][3]
double UpperCorner [2][3]
double Xmintot [2][3]
double Xmaxtot [2][3]
double TotalMeshSize [2]
double TimeLimitCPU
double CPU_TreeConstruction
double CPU_TreeWalk
double CPU_Gravity
double CPU_Potential
double CPU_Domain
double CPU_Snapshot
double CPU_Total
double CPU_CommSum
double CPU_Imbalance
double CPU_HydCompWalk
double CPU_HydCommSumm
double CPU_HydImbalance
double CPU_Hydro
double CPU_EnsureNgb
double CPU_Predict
double CPU_TimeLine
double CPU_PM
double CPU_Peano
double ErrTolTheta
double ErrTolForceAcc
double ErrTolIntAccuracy
double MinSizeTimestep
double MaxSizeTimestep
double MaxRMSDisplacementFac
double CourantFac
double TreeDomainUpdateFrequency
double MinGasHsmlFractional
double MinGasHsml
double SofteningGas
double SofteningHalo
double SofteningDisk
double SofteningBulge
double SofteningStars
double SofteningBndry
double SofteningGasMaxPhys
double SofteningHaloMaxPhys
double SofteningDiskMaxPhys
double SofteningBulgeMaxPhys
double SofteningStarsMaxPhys
double SofteningBndryMaxPhys
double SofteningTable [6]
double ForceSoftening [6]
double MassTable [6]
char InitCondFile [100]
char OutputDir [100]
char SnapshotFileBase [100]
char EnergyFile [100]
char CpuFile [100]
char InfoFile [100]
char TimingsFile [100]
char RestartFile [100]
char ResubmitCommand [100]
char OutputListFilename [100]
double OutputListTimes [500]
int OutputListLength


Detailed Description

This structure contains data which is the SAME for all tasks (mostly code parameters read from the parameter file). Holding this data in a structure is convenient for writing/reading the restart file, and it allows the introduction of new global variables in a simple way. The only thing to do is to introduce them into this structure.

Definition at line 246 of file allvars.h.


Field Documentation

double global_data_all_processes::ArtBulkViscConst
 

Sets the parameter $\alpha$ of the artificial viscosity

Definition at line 283 of file allvars.h.

Referenced by begrun(), hydro_evaluate(), and read_parameter_file().

double global_data_all_processes::Asmth[2]
 

Gives the scale of the long-range/short-range split (in mesh-cells), both for the coarse and the high-res mesh

Definition at line 368 of file allvars.h.

Referenced by pm_init_periodic(), pmforce_periodic(), and pmpotential_periodic().

double global_data_all_processes::BoxSize
 

Boxsize in case periodic boundary conditions are used

Definition at line 254 of file allvars.h.

Referenced by density(), do_box_wrapping(), fill_write_buffer(), force_treeevaluate_ewald_correction(), hydro_force(), ngb_treeallocate(), pm_init_periodic(), pmforce_periodic(), pmpotential_periodic(), read_parameter_file(), and write_file().

int global_data_all_processes::BufferSize
 

size of communication buffer in MB

Definition at line 264 of file allvars.h.

Referenced by allocate_commbuffers(), begrun(), read_file(), read_parameter_file(), and write_file().

int global_data_all_processes::BunchSizeDensity
 

number of particles fitting into the communication buffer in the density computation

Definition at line 266 of file allvars.h.

Referenced by allocate_commbuffers(), begrun(), and density().

int global_data_all_processes::BunchSizeDomain
 

number of particles fitting into the communication buffer in the domain decomposition

Definition at line 268 of file allvars.h.

Referenced by allocate_commbuffers(), begrun(), and domain_findExchangeNumbers().

int global_data_all_processes::BunchSizeForce
 

number of particles fitting into the buffer in the parallel tree-force algorithm

Definition at line 265 of file allvars.h.

Referenced by allocate_commbuffers(), begrun(), compute_potential(), gravity_forcetest(), and gravity_tree().

int global_data_all_processes::BunchSizeHydro
 

number of particles fitting into the communication buffer in the SPH hydrodynamical force computation

Definition at line 267 of file allvars.h.

Referenced by allocate_commbuffers(), begrun(), and hydro_force().

int global_data_all_processes::ComovingIntegrationOn
 

flags that comoving integration is enabled

Definition at line 323 of file allvars.h.

Referenced by advance_and_find_timesteps(), begrun(), compute_global_quantities_of_system(), compute_potential(), every_timestep_stuff(), fill_write_buffer(), find_next_outputtime(), gravity_forcetest(), gravity_tree(), hydro_evaluate(), hydro_force(), init(), long_range_force(), move_particles(), read_parameter_file(), readjust_timebase(), set_softenings(), and write_file().

double global_data_all_processes::Corner[2][3]
 

lower left corner of coarse and high-res PM-mesh

Definition at line 370 of file allvars.h.

double global_data_all_processes::CourantFac
 

SPH-Courant factor

Definition at line 422 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::CPU_CommSum
 

accumulated time used for communication, and for collecting partial results, in tree-gravity

Definition at line 388 of file allvars.h.

Referenced by every_timestep_stuff(), gravity_tree(), and main().

double global_data_all_processes::CPU_Domain
 

cumulative time spent for domain decomposition

Definition at line 385 of file allvars.h.

Referenced by domain_Decomposition(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_EnsureNgb
 

time needed to iterate on correct neighbour numbers

Definition at line 394 of file allvars.h.

Referenced by density(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_Gravity
 

cumulative time used for gravity computation (tree-algorithm only)

Definition at line 383 of file allvars.h.

Referenced by compute_accelerations(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_HydCommSumm
 

cumulative time used for communication in SPH, and for collecting partial results

Definition at line 391 of file allvars.h.

Referenced by density(), every_timestep_stuff(), hydro_force(), and main().

double global_data_all_processes::CPU_HydCompWalk
 

time used for actual SPH computations, including neighbour search

Definition at line 390 of file allvars.h.

Referenced by density(), every_timestep_stuff(), hydro_force(), and main().

double global_data_all_processes::CPU_HydImbalance
 

cumulative time lost due to work-load imbalance in SPH

Definition at line 392 of file allvars.h.

Referenced by density(), every_timestep_stuff(), hydro_force(), and main().

double global_data_all_processes::CPU_Hydro
 

cumulative time spent for SPH related computations

Definition at line 393 of file allvars.h.

Referenced by compute_accelerations(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_Imbalance
 

cumulative time lost accross all processors as work-load imbalance in gravitational tree

Definition at line 389 of file allvars.h.

Referenced by every_timestep_stuff(), gravity_tree(), and main().

double global_data_all_processes::CPU_Peano
 

time required to establish Peano-Hilbert order

Definition at line 398 of file allvars.h.

Referenced by domain_Decomposition(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_PM
 

time used for long-range gravitational force

Definition at line 397 of file allvars.h.

Referenced by compute_accelerations(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_Potential
 

time used for computing gravitational potentials

Definition at line 384 of file allvars.h.

Referenced by compute_potential(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_Predict
 

cumulative time to drift the system forward in time, including dynamic tree updates

Definition at line 395 of file allvars.h.

Referenced by compute_accelerations(), every_timestep_stuff(), main(), and move_particles().

double global_data_all_processes::CPU_Snapshot
 

time used for writing snapshot files

Definition at line 386 of file allvars.h.

Referenced by every_timestep_stuff(), main(), and savepositions().

double global_data_all_processes::CPU_TimeLine
 

time used for determining new timesteps, and for organizing the timestepping, including kicks of active particles

Definition at line 396 of file allvars.h.

Referenced by advance_and_find_timesteps(), every_timestep_stuff(), and main().

double global_data_all_processes::CPU_Total
 

cumulative time spent for domain decomposition

Definition at line 387 of file allvars.h.

Referenced by every_timestep_stuff(), main(), and run().

double global_data_all_processes::CPU_TreeConstruction
 

time spent for constructing the gravitational tree

Definition at line 381 of file allvars.h.

Referenced by compute_potential(), every_timestep_stuff(), gravity_tree(), and main().

double global_data_all_processes::CPU_TreeWalk
 

actual time spent for pure tree-walks

Definition at line 382 of file allvars.h.

Referenced by every_timestep_stuff(), gravity_tree(), and main().

char global_data_all_processes::CpuFile[100]
 

name of file with cpu-time statistics

Definition at line 468 of file allvars.h.

Referenced by begrun(), open_outputfiles(), and read_parameter_file().

double global_data_all_processes::CpuTimeBetRestartFile
 

cpu-time between regularly generated restart files

Definition at line 336 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and run().

double global_data_all_processes::DesNumNgb
 

Desired number of SPH neighbours

Definition at line 280 of file allvars.h.

Referenced by density(), read_parameter_file(), and setup_smoothinglengths().

char global_data_all_processes::EnergyFile[100]
 

name of file with energy statistics

Definition at line 467 of file allvars.h.

Referenced by begrun(), open_outputfiles(), and read_parameter_file().

double global_data_all_processes::ErrTolForceAcc
 

parameter for relative opening criterion in tree walk

Definition at line 403 of file allvars.h.

Referenced by begrun(), force_treeevaluate(), force_treeevaluate_potential(), force_treeevaluate_potential_shortrange(), force_treeevaluate_shortrange(), and read_parameter_file().

double global_data_all_processes::ErrTolIntAccuracy
 

accuracy tolerance parameter $ \eta $ for timestep criterion. The timestep is $ \Delta t = \sqrt{\frac{2 \eta eps}{a}} $

Definition at line 408 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::ErrTolTheta
 

BH tree opening angle

Definition at line 402 of file allvars.h.

Referenced by force_treeevaluate_ewald_correction(), gravity_tree(), and read_parameter_file().

double global_data_all_processes::ForceSoftening[6]
 

the same, but multiplied by a factor 2.8 - at that scale the force is Newtonian

Definition at line 453 of file allvars.h.

Referenced by force_treebuild_single(), force_treeupdate_pseudos(), force_update_node_recursive(), and set_softenings().

double global_data_all_processes::G
 

Gravity-constant in internal units

Definition at line 297 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_potential(), gravity_forcetest(), gravity_tree(), pmforce_periodic(), pmpotential_periodic(), and set_units().

double global_data_all_processes::GravityConstantInternal
 

If set to zero in the parameterfile, the internal value of the gravitational constant is set to the Newtonian value based on the system of units specified. Otherwise the value provided is taken as internal gravity constant G.

Definition at line 307 of file allvars.h.

Referenced by read_parameter_file(), and set_units().

double global_data_all_processes::Hubble
 

Hubble-constant in internal units

Definition at line 314 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_potential(), drift_integ(), gravity_forcetest(), gravity_tree(), gravkick_integ(), hydro_force(), hydrokick_integ(), init_drift_table(), long_range_force(), and set_units().

double global_data_all_processes::HubbleParam
 

little `h', i.e. Hubble constant in units of 100 km/s/Mpc. Only needed to get absolute physical values for cooling physics

Definition at line 318 of file allvars.h.

Referenced by read_parameter_file(), and write_file().

int global_data_all_processes::ICFormat
 

selects different versions of IC file-format

Definition at line 256 of file allvars.h.

Referenced by find_files(), init(), read_file(), read_ic(), and read_parameter_file().

char global_data_all_processes::InfoFile[100]
 

name of log-file with a list of the timesteps taken

Definition at line 469 of file allvars.h.

Referenced by begrun(), open_outputfiles(), and read_parameter_file().

char global_data_all_processes::InitCondFile[100]
 

filename of initial conditions

Definition at line 464 of file allvars.h.

Referenced by init(), and read_parameter_file().

double global_data_all_processes::InitGasTemp
 

may be used to set the temperature in the IC's

Definition at line 284 of file allvars.h.

Referenced by read_ic(), and read_parameter_file().

double global_data_all_processes::MassTable[6]
 

Table with particle masses for particle types with equal mass. If particle masses are all equal for one type, the corresponding entry in MassTable is set to this value, allowing the size of the snapshot files to be reduced.

Definition at line 456 of file allvars.h.

Referenced by get_particles_in_block(), read_file(), read_ic(), and write_file().

double global_data_all_processes::MaxNumNgbDeviation
 

Maximum allowed deviation neighbour number

Definition at line 281 of file allvars.h.

Referenced by begrun(), density(), and read_parameter_file().

int global_data_all_processes::MaxPart
 

This gives the maxmimum number of particles that can be stored on one processor.

Definition at line 251 of file allvars.h.

Referenced by allocate_memory(), density_evaluate(), domain_Decomposition(), domain_findExchangeNumbers(), force_treebuild_single(), force_treeevaluate_ewald_correction(), force_update_node_recursive(), free_memory(), gravity_tree(), hydro_evaluate(), init(), move_particles(), ngb_treefind_pairs(), ngb_treefind_variable(), pmforce_periodic(), pmpotential_periodic(), read_file(), and restart().

int global_data_all_processes::MaxPartSph
 

This gives the maxmimum number of SPH particles that can be stored on one processor.

Definition at line 252 of file allvars.h.

Referenced by allocate_memory(), domain_Decomposition(), domain_findExchangeNumbers(), free_memory(), read_file(), and restart().

double global_data_all_processes::MaxRMSDisplacementFac
 

this determines a global timestep criterion for cosmological simulations in comoving coordinates. To this end, the code computes the rms velocity of all particles, and limits the timestep such that the rms displacement is a fraction of the mean particle separation (determined from the particle mass and the cosmological parameters). This parameter specifies this fraction.

Definition at line 415 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::MaxSizeTimestep
 

maximum allowed timestep

Definition at line 413 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::MinEgySpec
 

the minimum allowed temperature expressed as energy per unit mass

Definition at line 286 of file allvars.h.

Referenced by advance_and_find_timesteps(), fill_write_buffer(), read_ic(), and set_units().

double global_data_all_processes::MinGasHsml
 

minimum allowed SPH smoothing length

Definition at line 435 of file allvars.h.

Referenced by density(), move_particles(), and set_softenings().

double global_data_all_processes::MinGasHsmlFractional
 

minimum allowed SPH smoothing length in units of SPH gravitational softening length

Definition at line 434 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::MinGasTemp
 

may be used to set a floor for the gas temperature

Definition at line 285 of file allvars.h.

Referenced by read_parameter_file(), and set_units().

double global_data_all_processes::MinSizeTimestep
 

minimum allowed timestep. Normally, the simulation terminates if the timestep determined by the timestep criteria falls below this limit.

Definition at line 411 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

int global_data_all_processes::NumCurrentTiStep
 

counts the number of system steps taken up to this point

Definition at line 340 of file allvars.h.

Referenced by every_timestep_stuff(), gravity_tree(), init(), and run().

int global_data_all_processes::NumFilesPerSnapshot
 

number of files in multi-file snapshot dumps

Definition at line 260 of file allvars.h.

Referenced by begrun(), read_parameter_file(), savepositions(), and write_file().

int global_data_all_processes::NumFilesWrittenInParallel
 

maximum number of files that may be written simultaneously when writing/reading restart-files, or when writing snapshot files

Definition at line 261 of file allvars.h.

Referenced by begrun(), read_ic(), read_parameter_file(), restart(), and savepositions().

long long global_data_all_processes::NumForcesSinceLastDomainDecomp
 

count particle updates since last domain decomposition

Definition at line 292 of file allvars.h.

Referenced by advance_and_find_timesteps(), domain_Decomposition(), init(), move_particles(), pmforce_periodic(), pmpotential_periodic(), and savepositions().

double global_data_all_processes::Omega0
 

matter density in units of the critical density (at z=0)

Definition at line 315 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_potential(), drift_integ(), gravity_forcetest(), gravity_tree(), gravkick_integ(), growthfactor_integ(), hydro_force(), hydrokick_integ(), long_range_force(), read_parameter_file(), and write_file().

double global_data_all_processes::OmegaBaryon
 

baryon density in units of the critical density (at z=0)

Definition at line 317 of file allvars.h.

Referenced by read_parameter_file().

double global_data_all_processes::OmegaLambda
 

vaccum energy density relative to crictical density (at z=0)

Definition at line 316 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_potential(), drift_integ(), gravity_forcetest(), gravity_tree(), gravkick_integ(), growthfactor_integ(), hydro_force(), hydrokick_integ(), long_range_force(), read_parameter_file(), and write_file().

char global_data_all_processes::OutputDir[100]
 

output directory of the code

Definition at line 465 of file allvars.h.

Referenced by begrun(), gravity_forcetest(), open_outputfiles(), read_parameter_file(), restart(), run(), and savepositions().

char global_data_all_processes::OutputListFilename[100]
 

name of file with list of desired output times

Definition at line 473 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

int global_data_all_processes::OutputListLength
 

number of output times stored in the table of desired output times

Definition at line 476 of file allvars.h.

Referenced by begrun(), find_next_outputtime(), read_outputlist(), and read_parameter_file().

int global_data_all_processes::OutputListOn
 

flags that output times are listed in a specified file

Definition at line 328 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::OutputListTimes[500]
 

table with desired output times

Definition at line 475 of file allvars.h.

Referenced by begrun(), find_next_outputtime(), and read_outputlist().

double global_data_all_processes::PartAllocFactor
 

in order to maintain work-load balance, the particle load will usually NOT be balanced. Each processor allocates memory for PartAllocFactor times the average number of particles to allow for that

Definition at line 270 of file allvars.h.

Referenced by read_file(), read_parameter_file(), and restart().

int global_data_all_processes::PeriodicBoundariesOn
 

flags that periodic boundaries are enabled

Definition at line 324 of file allvars.h.

Referenced by compute_potential(), init(), and read_parameter_file().

int global_data_all_processes::PM_Ti_begstep
 

end of present long-range timestep

Definition at line 361 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_global_quantities_of_system(), fill_write_buffer(), init(), and readjust_timebase().

int global_data_all_processes::PM_Ti_endstep
 

begin of present long-range timestep

Definition at line 360 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_accelerations(), compute_global_quantities_of_system(), domain_Decomposition(), fill_write_buffer(), gravity_forcetest(), init(), and readjust_timebase().

double global_data_all_processes::Rcut[2]
 

Gives the maximum radius for which the short-range force is evaluated with the tree (in mesh-cells), both for the coarse and the high-res mesh

Definition at line 369 of file allvars.h.

Referenced by pm_init_periodic().

char global_data_all_processes::RestartFile[100]
 

basename of restart-files

Definition at line 471 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and restart().

char global_data_all_processes::ResubmitCommand[100]
 

name of script-file that will be executed for automatic restart

Definition at line 472 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and run().

int global_data_all_processes::ResubmitOn
 

flags that automatic resubmission of job to queue system is enabled

Definition at line 325 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and run().

int global_data_all_processes::SnapFormat
 

selects different versions of snapshot file-formats

Definition at line 258 of file allvars.h.

Referenced by begrun(), read_parameter_file(), savepositions(), and write_file().

char global_data_all_processes::SnapshotFileBase[100]
 

basename to construct the names of snapshotf files

Definition at line 466 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and savepositions().

int global_data_all_processes::SnapshotFileCount
 

number of snapshot that is written next

Definition at line 333 of file allvars.h.

Referenced by init(), and run().

double global_data_all_processes::SofteningBndry
 

comoving gravitational softening lengths for type 5

Definition at line 443 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningBndryMaxPhys
 

maximum physical softening length for type 5

Definition at line 450 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningBulge
 

comoving gravitational softening lengths for type 3

Definition at line 441 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningBulgeMaxPhys
 

maximum physical softening length for type 3

Definition at line 448 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningDisk
 

comoving gravitational softening lengths for type 2

Definition at line 440 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningDiskMaxPhys
 

maximum physical softening length for type 2

Definition at line 447 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningGas
 

comoving gravitational softening lengths for type 0

Definition at line 438 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningGasMaxPhys
 

maximum physical softening length for type 0

Definition at line 445 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningHalo
 

comoving gravitational softening lengths for type 1

Definition at line 439 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningHaloMaxPhys
 

maximum physical softening length for type 1

Definition at line 446 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningStars
 

comoving gravitational softening lengths for type 4

Definition at line 442 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningStarsMaxPhys
 

maximum physical softening length for type 4

Definition at line 449 of file allvars.h.

Referenced by read_parameter_file(), and set_softenings().

double global_data_all_processes::SofteningTable[6]
 

current (comoving) gravitational softening lengths for each particle type

Definition at line 452 of file allvars.h.

Referenced by compute_potential(), domain_decompose(), and set_softenings().

int global_data_all_processes::Ti_Current
 

current time on integer timeline

Definition at line 354 of file allvars.h.

Referenced by advance_and_find_timesteps(), begrun(), compute_accelerations(), compute_global_quantities_of_system(), density(), domain_Decomposition(), fill_write_buffer(), gravity_forcetest(), gravity_tree(), hydro_force(), init(), move_particles(), readjust_timebase(), and run().

int global_data_all_processes::Ti_nextoutput
 

next output time on integer timeline

Definition at line 355 of file allvars.h.

Referenced by begrun().

double global_data_all_processes::Time
 

current time of the simulation

Definition at line 345 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_global_quantities_of_system(), energy_statistics(), every_timestep_stuff(), fill_write_buffer(), force_treebuild(), gravity_forcetest(), gravity_tree(), hydro_force(), init(), read_file(), restart(), run(), set_softenings(), and write_file().

double global_data_all_processes::Timebase_interval
 

factor to convert from floating point time interval to integer timeline

Definition at line 353 of file allvars.h.

Referenced by advance_and_find_timesteps(), compute_global_quantities_of_system(), density(), fill_write_buffer(), find_next_outputtime(), get_drift_factor(), get_gravkick_factor(), get_hydrokick_factor(), hydro_evaluate(), init(), move_particles(), and readjust_timebase().

double global_data_all_processes::TimeBegin
 

time of initial conditions of the simulation

Definition at line 346 of file allvars.h.

Referenced by find_next_outputtime(), init(), init_drift_table(), read_file(), read_parameter_file(), and readjust_timebase().

double global_data_all_processes::TimeBetSnapshot
 

simulation time interval between snapshot files

Definition at line 334 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::TimeBetStatistics
 

simulation time interval between computations of energy statistics

Definition at line 338 of file allvars.h.

Referenced by begrun(), init(), read_parameter_file(), and run().

double global_data_all_processes::TimeLastRestartFile
 

cpu-time when last restart-file was written

Definition at line 337 of file allvars.h.

Referenced by begrun(), and run().

double global_data_all_processes::TimeLastStatistics
 

simulation time when the energy statistics was computed the last time

Definition at line 339 of file allvars.h.

Referenced by init(), and run().

double global_data_all_processes::TimeLimitCPU
 

CPU time limit as defined in parameterfile

Definition at line 380 of file allvars.h.

Referenced by begrun(), read_parameter_file(), and run().

double global_data_all_processes::TimeMax
 

marks the point of time until the simulation is to be evolved

Definition at line 348 of file allvars.h.

Referenced by begrun(), find_next_outputtime(), init(), init_drift_table(), read_parameter_file(), readjust_timebase(), and run().

double global_data_all_processes::TimeOfFirstSnapshot
 

simulation time of first snapshot files

Definition at line 335 of file allvars.h.

Referenced by read_parameter_file().

double global_data_all_processes::TimeStep
 

difference between current times of previous and current timestep

Definition at line 347 of file allvars.h.

Referenced by every_timestep_stuff(), and gravity_tree().

char global_data_all_processes::TimingsFile[100]
 

name of file with performance metrics of gravitational tree algorithm

Definition at line 470 of file allvars.h.

Referenced by begrun(), open_outputfiles(), and read_parameter_file().

double global_data_all_processes::TotalMeshSize[2]
 

total extension of coarse and high-res PM-mesh

Definition at line 374 of file allvars.h.

long long global_data_all_processes::TotN_gas
 

total gas particle number (global value)

Definition at line 249 of file allvars.h.

Referenced by compute_accelerations(), domain_decompose(), read_file(), and restart().

long long global_data_all_processes::TotNumOfForces
 

counts total number of force computations

Definition at line 291 of file allvars.h.

Referenced by gravity_tree(), and init().

long long global_data_all_processes::TotNumPart
 

total particle numbers (global value)

Definition at line 248 of file allvars.h.

Referenced by advance_and_find_timesteps(), domain_Decomposition(), gravity_forcetest(), gravity_tree(), init(), move_particles(), pmforce_periodic(), pmpotential_periodic(), read_file(), read_ic(), restart(), and savepositions().

double global_data_all_processes::TreeAllocFactor
 

Each processor allocates a number of nodes which is TreeAllocFactor times the maximum(!) number of particles. Note: A typical local tree for N particles needs usually about ~0.65*N nodes.

Definition at line 274 of file allvars.h.

Referenced by gravity_tree(), init(), pmforce_periodic(), pmpotential_periodic(), read_parameter_file(), and restart().

double global_data_all_processes::TreeDomainUpdateFrequency
 

controls frequency of domain decompositions

Definition at line 427 of file allvars.h.

Referenced by advance_and_find_timesteps(), begrun(), domain_Decomposition(), init(), move_particles(), pmforce_periodic(), pmpotential_periodic(), read_parameter_file(), and savepositions().

int global_data_all_processes::TypeOfOpeningCriterion
 

determines tree cell-opening criterion: 0 for Barnes-Hut, 1 for relative criterion

Definition at line 326 of file allvars.h.

Referenced by begrun(), compute_accelerations(), gravity_tree(), and read_parameter_file().

int global_data_all_processes::TypeOfTimestepCriterion
 

gives type of timestep criterion (only 0 supported right now - unlike gadget-1.1)

Definition at line 327 of file allvars.h.

Referenced by begrun(), and read_parameter_file().

double global_data_all_processes::UnitCoolingRate_in_cgs
 

factor to convert internal cooling rate to cgs units

Definition at line 304 of file allvars.h.

Referenced by set_units().

double global_data_all_processes::UnitDensity_in_cgs
 

factor to convert internal length unit to g/cm^3*h^2

Definition at line 303 of file allvars.h.

Referenced by set_units().

double global_data_all_processes::UnitEnergy_in_cgs
 

factor to convert internal energy to cgs units

Definition at line 305 of file allvars.h.

Referenced by read_ic(), and set_units().

double global_data_all_processes::UnitLength_in_cm
 

factor to convert internal length unit to cm/h

Definition at line 301 of file allvars.h.

Referenced by read_parameter_file(), and set_units().

double global_data_all_processes::UnitMass_in_g
 

factor to convert internal mass unit to grams/h

Definition at line 299 of file allvars.h.

Referenced by read_ic(), read_parameter_file(), and set_units().

double global_data_all_processes::UnitPressure_in_cgs
 

factor to convert internal pressure unit to cgs units (little 'h' still around!)

Definition at line 302 of file allvars.h.

Referenced by set_units().

double global_data_all_processes::UnitTime_in_Megayears
 

factor to convert internal time to megayears/h

Definition at line 306 of file allvars.h.

Referenced by set_units().

double global_data_all_processes::UnitTime_in_s
 

factor to convert internal time unit to seconds/h

Definition at line 298 of file allvars.h.

Referenced by set_units().

double global_data_all_processes::UnitVelocity_in_cm_per_s
 

factor to convert intqernal velocity unit to cm/sec

Definition at line 300 of file allvars.h.

Referenced by read_parameter_file(), and set_units().

double global_data_all_processes::UpperCorner[2][3]
 

upper right corner of coarse and high-res PM-mesh

Definition at line 371 of file allvars.h.

double global_data_all_processes::Xmaxtot[2][3]
 

maximum particle coordinates both for coarse and high-res PM-mesh

Definition at line 373 of file allvars.h.

double global_data_all_processes::Xmintot[2][3]
 

minimum particle coordinates both for coarse and high-res PM-mesh

Definition at line 372 of file allvars.h.


The documentation for this struct was generated from the following file:
Generated on Sun May 22 17:33:30 2005 for GADGET-2 by  doxygen 1.3.9.1