GADGET-4
simulation.h
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#ifndef SIMULATION_H
13#define SIMULATION_H
14
15#include "gadgetconfig.h"
16
17#include <mpi.h>
18
19#include "../cooling_sfr/cooling.h"
20#include "../data/allvars.h"
21#include "../data/constants.h"
22#include "../data/dtypes.h"
23#include "../data/lcparticles.h"
24#include "../data/macros.h"
25#include "../data/mmparticles.h"
26#include "../data/mymalloc.h"
27#include "../data/simparticles.h"
28#include "../domain/domain.h"
29#include "../fmm/fmm.h"
30#include "../fof/fof.h"
31#include "../gravity/ewald.h"
32#include "../gravity/grav_forcetest.h"
33#include "../gravtree/gravtree.h"
34#include "../gravtree/gwalk.h"
35#include "../io/parameters.h"
36#include "../io/restart.h"
37#include "../io/test_io_bandwidth.h"
38#include "../lightcone/lightcone.h"
39#include "../logs/logs.h"
40#include "../mergertree/mergertree.h"
41#include "../mpi_utils/mpi_utils.h"
42#include "../mpi_utils/setcomm.h"
43#include "../ngbtree/ngbtree.h"
44#include "../ngenic/ngenic.h"
45#include "../pm/pm.h"
46#include "../sph/sph.h"
47#include "../system/pinning.h"
48
49class sim : public pinning, public test_io_bandwidth
50{
51 public:
52 sim(MPI_Comm comm) : setcomm(comm), test_io_bandwidth(comm) {}
53
54 /* here come the main classes the code operates on */
55
56 simparticles Sp{Communicator}; /* stores all the simulation particles of the simulation */
57
58 domain<simparticles> Domain{Communicator, &Sp}; /* get an instance of a domain decomposition, operating on Sp */
59
60 /* Note: GravTree/NgbTree inherit their communicator and NTask/ThisTask from their associated domain object */
61
62#ifdef FMM
63 fmm GravTree; /* get an instance of a gravitational tree */
64#else
65 gwalk GravTree; /* get an instance of a gravitational tree */
66#endif
67
68 sph NgbTree; /* get an instance of a neighbour search tree */
69
70#ifdef PMGRID
71 pm PM{Communicator};
72#endif
73
74#ifdef MERGERTREE
75 mergertree MergerTree{Communicator, &Sp};
76#endif
77
78#ifdef COOLING
79 coolsfr CoolSfr{Communicator};
80#endif
81
82#ifdef NGENIC
83 ngenic Ngenic{Communicator, &Sp};
84#endif
85
86#ifdef LIGHTCONE
87#ifdef LIGHTCONE_PARTICLES
88 lcparticles Lp{Communicator}; /* stores all the buffered light-cone particles of the simulation */
89#endif
90#ifdef LIGHTCONE_MASSMAPS
91 mmparticles Mp{Communicator}; /* stores buffered particles to construct projected mass maps on the lightcone */
92#endif
93
94#if defined(LIGHTCONE_PARTICLES) && defined(LIGHTCONE_MASSMAPS) /* both particles and massmaps */
95 lightcone LightCone{Communicator, &Sp, &Lp, &Mp};
96#else
97#if defined(LIGHTCONE_PARTICLES)
98 lightcone LightCone{Communicator, &Sp, &Lp};
99#else
100 lightcone LightCone{Communicator, &Sp, &Mp};
101#endif
102#endif
103#endif // end of LIGHTCONE
104
105#ifdef FORCETEST
106 gravtest GravTest{&Sp, &GravTree, &Domain};
107#endif
108
109 void rearrange_lightcone(int argc, char **argv);
110 void rearrange_snapshot(int argc, char **argv);
111
112 template <typename partset>
113 void rearrange_generic(partset &Tp, int conenr, int firstnr, int lastnr);
114
115 template <typename partset>
116 void rearrange_fill_treetable(partset &Tp);
117
118 template <typename partset>
119 void rearrange_read(partset &Tp, int num, int conenr);
120
121 template <typename partset>
122 void rearrange_write(partset &Tp, int num, int conenr);
123
124 private:
125#ifdef PERIODIC
126 void check_omega(void);
127#endif
128
129 void setup_smoothinglengths(void);
130 void recreate_unique_ids(void);
131 void test_id_uniqueness(void);
132 int check_for_interruption_of_run(void);
133 void set_non_standard_physics_for_current_time(void);
134 void calculate_non_standard_physics_end_of_step(void);
135 integertime find_next_outputtime(integertime ti_curr);
136 double measure_cpu_performance(MPI_Comm Communicator);
137 double measure_hyper_cube_speed(const char *tag, MPI_Comm Communicator);
138 void measure_iprobe_performance(const char *tag);
139
140#ifdef SECOND_ORDER_LPT_ICS
141 double F1_Omega(double a);
142 double F2_Omega(double a);
143 int executed = 0;
144 void second_order_ic_correction(void);
145#endif
146
147 public:
148 void hello(void);
149 void endrun(void);
150 void begrun1(const char *parameterFile);
151 void begrun2(void);
152 void init(int RestartSnapNum);
153 void run(void);
154 void set_units(void);
156 void healthtest(void);
158 long long report_comittable_memory(long long *MemTotal, long long *Committed_AS, long long *SwapTotal, long long *SwapFree);
159 long long report_free_size_in_tmpfs(void);
162 void do_hydro_step_first_half(void);
163 void do_hydro_step_second_half(void);
165 void find_hydro_timesteps(void);
166 void gravity(int timebin);
168 void gravity_comoving_factors(int timebin);
169 void gravity_pm(int timebin);
170 void gravity_set_oldacc(int timebin);
171
172 void hydro_force(int step_indicator);
173 void compute_grav_accelerations(int timebin);
174#ifdef FORCETEST_TESTFORCELAW
175 void gravity_forcetest_testforcelaw(void);
176#endif
177
178#ifdef EXTERNALGRAVITY
179 void gravity_external(void);
180#endif
181};
182
183#endif
Definition: domain.h:31
Definition: gwalk.h:18
MPI_Comm Communicator
Definition: setcomm.h:31
Definition: simulation.h:50
void rearrange_write(partset &Tp, int num, int conenr)
void run(void)
Definition: run.cc:49
void compute_grav_accelerations(int timebin)
This routine computes the gravitational accelerations for all active particles.
Definition: gravity.cc:41
sim(MPI_Comm comm)
Definition: simulation.h:52
void rearrange_generic(partset &Tp, int conenr, int firstnr, int lastnr)
domain< simparticles > Domain
Definition: simulation.h:58
void do_hydro_step_second_half(void)
Definition: kicks.cc:469
long long report_free_size_in_tmpfs(void)
Definition: system.cc:203
void hydro_force(int step_indicator)
Definition: kicks.cc:480
void begrun1(const char *parameterFile)
This function performs the initial set-up of the simulation.
Definition: begrun.cc:130
void gravity_comoving_factors(int timebin)
Definition: gravity.cc:277
gwalk GravTree
Definition: simulation.h:65
void mpi_report_comittable_memory(void)
Definition: system.cc:215
long long report_comittable_memory(long long *MemTotal, long long *Committed_AS, long long *SwapTotal, long long *SwapFree)
Definition: system.cc:167
void gravity(int timebin)
main driver routine of gravity tree/fmm force calculation
Definition: gravity.cc:171
void gravity_pm(int timebin)
void do_hydro_step_first_half(void)
Definition: kicks.cc:449
void find_hydro_timesteps(void)
Definition: timestep.cc:39
void rearrange_snapshot(int argc, char **argv)
void hello(void)
Definition: begrun.cc:58
void init(int RestartSnapNum)
Prepares the loaded initial conditions for the run.
Definition: init.cc:51
void gravity_long_range_force(void)
void set_units(void)
Computes conversion factors between internal code units and the cgs-system.
Definition: begrun.cc:335
void rearrange_lightcone(int argc, char **argv)
void find_timesteps_and_do_gravity_step_first_half(void)
performs the first half step kick operator for the gravity
Definition: kicks.cc:40
void rearrange_read(partset &Tp, int num, int conenr)
sph NgbTree
Definition: simulation.h:68
void gravity_set_oldacc(int timebin)
Definition: gravity.cc:248
void rearrange_fill_treetable(partset &Tp)
void find_global_timesteps(void)
void do_gravity_step_second_half(void)
performs the second gravity half step kick operator
Definition: kicks.cc:274
void begrun2(void)
This function does late setup, after the IC file has been loaded but before run() is called.
Definition: begrun.cc:262
void healthtest(void)
Definition: healthtest.cc:36
void create_snapshot_if_desired(void)
Check if a snapshot should be saved.
Definition: run.cc:507
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
Definition: sph.h:21
int integertime
Definition: constants.h:331